{
 "cells": [
  {
   "cell_type": "markdown",
   "id": "dfbedf1d",
   "metadata": {},
   "source": [
    "# rf314_paramfitrange\n",
    "Multidimensional models: working with parameterized ranges in a fit.\n",
    "This an example of a fit with an acceptance that changes per-event\n",
    "\n",
    "`pdf = exp(-t/tau)` with `t[tmin,5]`\n",
    "\n",
    "where `t` and `tmin` are both observables in the dataset\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:31 PM.</small></i>"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 1,
   "id": "fb125a9d",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:31:16.796989Z",
     "iopub.status.busy": "2026-05-19T20:31:16.796855Z",
     "iopub.status.idle": "2026-05-19T20:31:17.746622Z",
     "shell.execute_reply": "2026-05-19T20:31:17.745993Z"
    }
   },
   "outputs": [],
   "source": [
    "import ROOT"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "9e1e9e1a",
   "metadata": {},
   "source": [
    "Define observables and decay pdf\n",
    "---------------------------------------------------------------"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "9526c5db",
   "metadata": {},
   "source": [
    "Declare observables"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 2,
   "id": "fa209daa",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:31:17.754182Z",
     "iopub.status.busy": "2026-05-19T20:31:17.754049Z",
     "iopub.status.idle": "2026-05-19T20:31:17.930795Z",
     "shell.execute_reply": "2026-05-19T20:31:17.928854Z"
    }
   },
   "outputs": [],
   "source": [
    "t = ROOT.RooRealVar(\"t\", \"t\", 0, 5)\n",
    "tmin = ROOT.RooRealVar(\"tmin\", \"tmin\", 0, 0, 5)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "82aef061",
   "metadata": {},
   "source": [
    "Make parameterized range in t : [tmin,5]"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 3,
   "id": "f2404ae4",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:31:17.932442Z",
     "iopub.status.busy": "2026-05-19T20:31:17.932321Z",
     "iopub.status.idle": "2026-05-19T20:31:18.054111Z",
     "shell.execute_reply": "2026-05-19T20:31:18.053576Z"
    }
   },
   "outputs": [],
   "source": [
    "t.setRange(tmin, ROOT.RooFit.RooConst(t.getMax()))"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "0f9149c6",
   "metadata": {},
   "source": [
    "Make pdf"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 4,
   "id": "99bbe3fd",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:31:18.058502Z",
     "iopub.status.busy": "2026-05-19T20:31:18.058359Z",
     "iopub.status.idle": "2026-05-19T20:31:18.180235Z",
     "shell.execute_reply": "2026-05-19T20:31:18.179556Z"
    }
   },
   "outputs": [],
   "source": [
    "tau = ROOT.RooRealVar(\"tau\", \"tau\", -1.54, -10, -0.1)\n",
    "model = ROOT.RooExponential(\"model\", \"model\", t, tau)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "120352c4",
   "metadata": {},
   "source": [
    "Create input data\n",
    "------------------------------------"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "cfc0ad99",
   "metadata": {},
   "source": [
    "Generate complete dataset without acceptance cuts (for reference)"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 5,
   "id": "16832ba3",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:31:18.182323Z",
     "iopub.status.busy": "2026-05-19T20:31:18.182163Z",
     "iopub.status.idle": "2026-05-19T20:31:18.333114Z",
     "shell.execute_reply": "2026-05-19T20:31:18.332760Z"
    }
   },
   "outputs": [],
   "source": [
    "dall = model.generate({t}, 10000)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "f372e6d8",
   "metadata": {},
   "source": [
    "Generate a (fake) prototype dataset for acceptance limit values"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 6,
   "id": "de59a4db",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:31:18.337652Z",
     "iopub.status.busy": "2026-05-19T20:31:18.337508Z",
     "iopub.status.idle": "2026-05-19T20:31:18.454805Z",
     "shell.execute_reply": "2026-05-19T20:31:18.454162Z"
    }
   },
   "outputs": [],
   "source": [
    "tmp = ROOT.RooGaussian(\"gmin\", \"gmin\", tmin, 0.0, 0.5).generate({tmin}, 5000)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "3ee3eae7",
   "metadata": {},
   "source": [
    "Generate dataset with t values that observe (t>tmin)"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 7,
   "id": "0d90d731",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:31:18.456926Z",
     "iopub.status.busy": "2026-05-19T20:31:18.456804Z",
     "iopub.status.idle": "2026-05-19T20:31:18.586345Z",
     "shell.execute_reply": "2026-05-19T20:31:18.585708Z"
    }
   },
   "outputs": [],
   "source": [
    "dacc = model.generate({t}, ProtoData=tmp)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "20030120",
   "metadata": {},
   "source": [
    "Fit pdf to data in acceptance region\n",
    "-----------------------------------------------------------------------"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 8,
   "id": "787f8b1c",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:31:18.588428Z",
     "iopub.status.busy": "2026-05-19T20:31:18.588305Z",
     "iopub.status.idle": "2026-05-19T20:31:18.773279Z",
     "shell.execute_reply": "2026-05-19T20:31:18.772647Z"
    }
   },
   "outputs": [
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "[#1] INFO:Fitting -- RooAbsPdf::fitTo(model_over_model_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 751.474 μs\n",
      "[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_model_over_model_Int[t]_modelData) Summation contains a RooNLLVar, using its error level\n",
      "[#1] INFO:Minimization -- [fitFCN] No discrete parameters, performing continuous minimization only\n"
     ]
    }
   ],
   "source": [
    "r = model.fitTo(dacc, Save=True, PrintLevel=-1)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "2a8e4bcd",
   "metadata": {},
   "source": [
    "Plot fitted pdf on full and accepted data\n",
    "---------------------------------------------------------------------------------"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "7469fb37",
   "metadata": {},
   "source": [
    "Make plot frame, datasets and overlay model"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 9,
   "id": "26466bf6",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:31:18.775338Z",
     "iopub.status.busy": "2026-05-19T20:31:18.775210Z",
     "iopub.status.idle": "2026-05-19T20:31:18.961823Z",
     "shell.execute_reply": "2026-05-19T20:31:18.961217Z"
    }
   },
   "outputs": [
    {
     "data": {
      "text/plain": [
       "<cppyy.gbl.RooPlot object at 0x55e3fb6ded30>"
      ]
     },
     "execution_count": 9,
     "metadata": {},
     "output_type": "execute_result"
    },
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "[#1] INFO:Plotting -- RooPlot::updateFitRangeNorm: New event count of 5000 will supersede previous event count of 10000 for normalization of PDF projections\n"
     ]
    }
   ],
   "source": [
    "frame = t.frame(Title=\"Fit to data with per-event acceptance\")\n",
    "dall.plotOn(frame, MarkerColor=\"r\", LineColor=\"r\")\n",
    "model.plotOn(frame)\n",
    "dacc.plotOn(frame)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "3b3f8071",
   "metadata": {},
   "source": [
    "Print fit results to demonstrate absence of bias"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 10,
   "id": "f6e636e6",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:31:18.963546Z",
     "iopub.status.busy": "2026-05-19T20:31:18.963426Z",
     "iopub.status.idle": "2026-05-19T20:31:19.205523Z",
     "shell.execute_reply": "2026-05-19T20:31:19.204439Z"
    }
   },
   "outputs": [
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "\n",
      "  RooFitResult: minimized FCN value: 2823.97, estimated distance to minimum: 3.17108e-08\n",
      "                covariance matrix quality: Full, accurate covariance matrix\n",
      "                Status : MINIMIZE=0 HESSE=0 \n",
      "\n",
      "    Floating Parameter  InitialValue    FinalValue +/-  Error     GblCorr.\n",
      "  --------------------  ------------  --------------------------  --------\n",
      "                   tau   -1.5400e+00   -1.5335e+00 +/-  2.22e-02  <none>\n",
      "\n"
     ]
    },
    {
     "name": "stderr",
     "output_type": "stream",
     "text": [
      "Info in <TCanvas::Print>: png file rf314_paramranges.png has been created\n"
     ]
    }
   ],
   "source": [
    "r.Print(\"v\")\n",
    "\n",
    "c = ROOT.TCanvas(\"rf314_paramranges\", \"rf314_paramranges\", 600, 600)\n",
    "ROOT.gPad.SetLeftMargin(0.15)\n",
    "frame.GetYaxis().SetTitleOffset(1.6)\n",
    "frame.Draw()\n",
    "\n",
    "c.SaveAs(\"rf314_paramranges.png\")"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "9ff83f6a",
   "metadata": {},
   "source": [
    "Draw all canvases "
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 11,
   "id": "ed13af2a",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:31:19.207075Z",
     "iopub.status.busy": "2026-05-19T20:31:19.206940Z",
     "iopub.status.idle": "2026-05-19T20:31:19.388752Z",
     "shell.execute_reply": "2026-05-19T20:31:19.388152Z"
    }
   },
   "outputs": [
    {
     "data": {
      "text/html": [
       "\n",
       "\n",
       "<div id=\"root_plot_1779222679379\" style=\"width: 600px; height: 600px; position: relative\">\n",
       "</div>\n",
       "\n",
       "</div>\n",
       "<script>\n",
       "   function process_root_plot_1779222679379() {\n",
       "      function execCode(Core) {\n",
       "         Core.settings.HandleKeys = false;\n",
       "         \n",
       "Core.unzipJSON(40753,'WkwIfjIAMZ8AeAHtneuTHDeO4P8VRcV+uIug8wi+kflJlq3z3Pmh8GMs7dyEoiRVS3VudWu7S7Y8G/7fL35gVnU9WpY8Y8969yxP9SSYJAG+ABAAmf++eLz58dXqYvlytRgXX99bXny/vP529eSri+Wr6xeXm4VbnH1zsf6316s/fbQYvVucfbjeXPenL57839XTDekLsn3xarO+vJiB/72+eLYYo1uc7Woa//02XG9DEGMKMSW3OPt0fbG6d3l+ebUYZQa/2vx4vroBv10/27zo4P31+fmcGWIBt5m9t+Krs81ny6vn64vF6AfJbnH25fr5i/00Gnm52Vy+3MvnFmdfX746THh4toaK4BZnj24e7/ZHSjy83iw3oFElzwF0t0MUvn+1fLk6Jpy0o5bv8h22aJd1m7ytkuL7XUOVH15ePVtdfbX+29x9e4mfXT5b9YF9KIvxAxl8/xe1xpIYykeW7mVIWnLUUKSGSjvDYsxD6f9qKd5Xkh+FxagSh1Ra1dBya3TDw83l3SfXD9ZvVuffLcamQwJNaK3EnOYM7/HWikblX4xRpdHsRwd1Z0lDJUPJ6oswlx5tdnX/zNvF+IEfcom5FJ+Kxlw9hb85qJ3p9c1NdfvgYsxa3OLszwcFcmWw/3xT5DBhMX7QE7a9s7l8+N3t42Dd8+7XzDwvIi2XFpJXLdLcYq/+R9/dNkBz7bz84JbB3r420mpNSZJKCZoqlT88v/zh84/u9Xn0aB94+M0re2HjtPf87S71k93T3SfXBxXdfXJ9UNfdJ9c3xe4+ub4p+c2bl7ayGZwfbx7fvFy+WYws9m9+tMcmMmgNIcesVQqz++sXq81yMUbKPnixnp/uXr9aPd18udysL3ujPn/98snqqj9/vX763Zubxx/746eXz+fETy+f36T9rb99sHz2YLm+gC+4xdm9q8vr6xfL9VzhDnxwOTPZ/QXL/OnwzVr97PLZ+my9erYYz5bn1yu3OPufV+tnbw7BH2/Au0+u711eXu3l//jZerN8AifbXL2mgvvrN6tnB+3eVv3gav1yvVl/v7o+4eafrq8RFltBMoPLq6vF+Je/usXlqw0PP7nF2cdvVk+vF+PF6/Nztzj7vEueq7Mo6fGr5dXy5dXy4vnK5MLX6w1k3f7y89cvHyzPV5vNVjjQnZ+v3mxOUz/601cPPr37aDEu/mX76BZnH12+fnK++vD12dl2OL9cbZbrC/py7omH1+u/rb653r5/dAja2y9Xy/PFGEBur/fgb9cXzy5/+PryFSvRLW7gR/vwzKFvMnyyQh7N0+OHLS+592IxM4x7y83mpPvvbjZddtOyhx+uNj+sVhezbDqArE/vX12+/Pry1WKUgTn18NlyA+s34NEWYE3c7YD85BbffXb5/eqLV8t/e72bKd99uaJHDhPPPlk/f/EpTZjlsM3a5ebpi223fvfVi8sfPv5+dbH5arPcvL7eTc7v7r7eXDI9djk/W128/nB51WEmz92nTL9dibMvV8tnX1yc/7gtcfbtevPi8vVmf6ZuZ+8ny+t57m1T9nP95Ug9+dWUIETHW5Wgb1dPjBusL56/TRNiZtw7X15fz0uFfF312k94xSRd+DHk7ObfJKN33vkpWCpPcfT21k/J3oWcp7wr46cy5yS9znl5bmNTF0RcS5Py3JILUiYBYXLzbxIZRYOTFpzUMEkYa3H2v0niKD64+TdJGiVUN/8myaPk6ObfJGWU2tz8m6SOwSc3/yZpYwjq5t8kOoYU3PybgrfMopT3U5BDMIyi0YmKk9KmEEdpxaE9SIlTSB1vyU6yn0IepUXeOMl5ClA1E5n8FOooOTtJwUn0U2ijxGwNbGUKOkrNLvgCJVP0o8TgDHWJU5RRYoUCJ8VPMYwSvZPcHBXH2PG01svSVaHjyXGKeRSpzsqnPMUyik9QYCTHOopXZwVimGIbJSTr1+B1inrT6zlPyY9SvbPmqkxJRobYaI55SmEMQZz1DyBd1TqR0U8pjdbFjFpIU8qj9SlkiE6pjFSjyWmZUu1o6KhUp9RGwyhG5pR0lFSdNULClL0RYVNsyjJKKs7m7JSD9cMMxFHb9jmNNW2f85hl+1zGULfPdYzbx8bEnqvR0bvEBM5T8aN3pfVnlo5Kfw4jjbJlNRVWkETtb1hEzA4rngFa7UABh3QkxVZSnLEYdlukMhXQGxDqVMFvgE9ThQAAaX6qUGBALlPdrmHGt0KBvfFtqlAA0PJUjYCcXZGpbldyrFPdopc4VR07hc5PzY+dQJ5lzNYMnsPYrBt4Zgl37uKnxqy0bgDoyxd8fmosk7IFWLqdKj81Fq51AwDL1nrb+UlnxhUbgHSOJAmgMy9R2JdGeyOlAKQOpADQuRiTyk/a2ZhazdVyVSG92XOqPKs9h+j8JN6PIuJSNqYgXmyGpuBSmMSHzhSU4ZVJfBxbdWKzhtfJ2iG5uAiYxxCai9XlNIkvo2TvKgtUJ/GwjGJ5pVAVyzO4GF0IgHSJK/SqTGK8NbOGwiQCGxMXxQWgYL0KMyo6CbyVZ34hTwJzZTTpeOqV3Cnu81ngrl6cqhPvJ5E6wuBoj7RJ4K7Mu6JO4Og6SlMQyyTBj5LgWclJgb/Dx3q/NcAwWjFfXdNJQhxDNDx+Elhry4gFxzIXeKtEZ6wskLnYUg/SOs1w1z0RELyMKixJlSkwJLDE/puChzdFN/+mwJgI7bHfFHzu/Ba2PgUGpUVjAX4Kvo6szSquUnMbQ9ROltQpeB0pBbPOU2BImAUUnoLIWMX1wlOQYHx5Lj0FBgV8FJ6CJEr2wlOQTMleeAqMBxxdKqWnILVTu0Xado0zrGr0btEi7mBHnegp2JBAbcfLkPS2SotTsCGxxoaYphDSrrWGmUGZWwvRAS5+06uhbtvbC7ddgykct2sY3ubTFLfLmMmE3GIlz1K8w30x23sfpsh63r43eE8nMLivavJnmSILe85uYNtpEAbC27q6kWVK/kYx8lPyXSegOFDYlQS60VOAEuIEiTIlT2907kX7ki1m4/3WnsRiBkYa5zglpo4Pu/Yl8R3mvQ9TQlmCRCSowQFcc/OSRIN8MvrF6CArrRETchQ2aMvq/ZTEWH1vl8DqKUI6goYCfkoBOdObFOYh6TrjlFCH5iYh0hKLdm6SNZn5MTepv4eTBGtyf4/2MTeJLght26SeXbdtstzRA4LAMkcxsI/alOJW6FnRaFKvj/eUogneXUGTvLtydIcN2pSZkvMcoZasshvqDt8MPRRlvRn8/j719vQ2Tpkp2fEYyVkLJPdRm7LWLWR1Wdt7A6es1nQgH6biGQKTilPxXUunTPG02QZtKp4W26BNxXdtHY2BWdhXmJFYYF99hc1wsybPwzYVZArv+7BNxViYzcCOklnYZ6A1qSBabIX16iRadTQRApEsNgnnzL0/ts2S/e1DmbV2uB4l++qYG4dUYTb27gg2ERi1qTAlt+2jHBxrbh+zpMCy+gqb4Vnwz8NfjGf19llXm5ZuK2yGu3jrE6QE3TYP3GjpvXUG9fU554w3m6gSe6egkRS0c9aYPffu6M82NXisW0bYyZ7qzAjnbp8qs24mMctUmXQzSOG8a6+9vBl+mldn9cZImOqs4IDJMvfW9ZdtywONOzR44NyvWabmOyewQZnazASNSTS/G3PezOR03tPm2TjTO7XtbDRO05iL89wzHLvORsfsHLmjmCchmcmIKJ3ZCjmNE85V3jBC3szsqDe3IUe3DNmHqSFH+7S3zG0HGY4dL0LJtQkIeoBZOvTemRoidMvXqRYmyWzrzLjBI2eQwrMmPJM0z8a5F6Y2z8Ytpr5WgSCJPeOMCGV75l2dKCbj3DaWRWPLOHcRy6LF3kl93kzthkFCcTR5AU0978y1Z0Hd4nZ/0DPf8GxDZFuUbdG0v0dpyTg29Vi9yTh2X3xTS7ZTmudJ2jJs2pluyIG6dENOr2ZHDg1LO3J65i05VJS31Bgw91AX0VNjuziPizUk76aR1ZuhaObOre8a5w6dWu4UzR0+tQxF9hZTSDaKbkAomvljsz3kLA9a6QKtD9vU+j6yD9rUbCM5M87WN5Kd2qmVLt77VJha2ROy4Lft5E2tvX+2rbYd5a5o31Ju39XZFLIlqM7Sb0tS7RN72ym1T+zOsVvtE3vb0NqVoV0vVHY22x12q31ibzuCjWZ/yYxqfWJ38ddaH7Zt1tZHrQuN1mzMupRozXSgmRx2nDa3+uCy55x5pvVR60t/i6Mv/R0OW/ozCp1ljxWbOfW2B2ZGPeOflQMTU21m0jPQmeJ2NGcevRvPmUdvq9kJIJvQM4+2inRm0TNww6HpN91j0cxghUnPstiHSWc23bHozKfnluiWT89DrTOf7p2ps86whbrGMI+0zox6fjlb12xOQRTbDTq+/yaFWTKP+2/SYCPYZaiGNObItloiRW2aJ7bGddJQxuxtB57Y9NcxKRt1SVTaxtQcG89UJg06puKyd5LapNGPCfbpJMukUcaUnBGTJo1hTNEZe6qTxjim4Nh+FywJaUzeYXErYdKYx6guK7awSWMZY3Olmwg11jFWxxa5gqGNMTs2wBUMitmkFCe1TJr8GKMr1UnVSZOMMbCbZ2+sKYzRtnzS8qQpjkEde+FWJ01pDM2MBIpZI2O2qhlr4aRY0oqr1YmCoWIhqOpEwcBm1TVmskyadAzBta7XqdnRnG396qQwR3UtY5OaFNbYXCsuSJw0Yz1wrbkgedKcRilOsV61SbNtCDGEBZkUvpidhhliyZt5IlAn+h0mlw7ZtllLhwp6ldM6Q0wZp+zzMdVgM8AWEUQnLRFIMFwYaMTYrDUw97dxfls62A0XWmoHaWSbtDSaJZ5WAmKQdZgmDKyYPh2m0iB10ioGYg4zEJOPwyTUQeshkTCDWIkd1gtM4ApnbE7M1gCIdaHbCw2sHQQvb9nMO8HKQW/DGQEhI0/avIFYCgyUDkIGb7Fdm9Gjg7GDkJEmxSDHW5oPiEHWSaD5gKWDkAFYOwgZcVIYJJkhA1ANjJARJ9VOFbYnCZOyiWpYH2ew9xX2JpFJ4ZJYpyEDsI9ghAzAPoIRMrDB2YSSCBmANqPMumqgTSlJzG+dlC19coIVy2NLw2gHzIRvwBijYBhbGKu6k4Q5wt7HUbBAFxd8JT8WTCzSWzh3uG1hzDJOkmIRJ381OEOOwUgV2I4LPvPezKmSoQdY/KizKDZQDDQ3Bm/DiDG920rESzQQWuxtMhBSIm/ziEkHSgwsgAVCsDhiwMOeuAWbgZBhb3XEmkevYGIM3sqadRBQQAQ7U/o0BIjE/aF0mam6Zm4VpcvQdcMM0ya0XWyQcCfINiPeDEOoqbvdzaHQgskIyuFfRkyf/sYjQc8GTHHbOGngj2gsMwx+NmG4lswrAYyFsdsUG/REbLEzDD2RzWG3i2Kf9NE8AmYn7ZbcCuMVzLON3ogNtozjBg4uHvYOfmDogcGDv/N78bB48AODHyYPfiQA+FNEWHQ7Jv2REqKkw9CTMoJG4OcV/KkghqTR5+BPFSHVYcPfEGHdW0V/JEXAdZj+yB7xZ+63Aj1ZEI5mFi3gzwFB2mHwzzqx+bgYj5zGTHkEOvjRi6uTxsYB/LmMyEiDDV8dC/Tw3vC1kbljMO3POuLiwu5rpvLisSB2GHxFRmzWvE/gK+Yw7DD4ShzpGzx0CXwljbgWDQZfySO4DAZfKWboNpj2lTqa1R49wvC1EecfsJnji46K9w8YfNXjVeow/YtY8BAMvwMhgoENjiWA0XRmusg4IglpRIz0BHAiHlg0Tc1NKB4BgZWVhGBYYSygJcHQ4hYALQmGFr8da9Sbs008giLgWiSBkUVUMHUtAbQIC/M4koO5hbhgMuL5xM7vERjMRnIIaBEZTL+bBGzyeC3JQf8jNphglgBhCI4004G3wCM6mGKWA8IQHmhnlkAHbW1wlgClZoSb+0Og1Owhc38YpzchQp/isYB0EyNznxovR5DYzCMHlCJKYEuWAKUIk85ncDWSgEIwjy0J5gQyTtCcwGbNDcTStgS8LltPkCWYhwa/w7YOXCnIFHNENAc/xRtksxkY34q5HmhaczBUcwjREa3iKBVzCJkDvToYLB4hUeZbdTBYfEKiNKs6GCxuIcGPA2xeIsw4ELSFu6XKYAg2NwRNQsUDZucCPcU8NriG8FqSHwaMbwiXsMHQh3PIJnPDySh4h3Cy0B4YLg6iYHO7uQZ92FLg+LhzwY91b57qMFxzEcHh1TsYLl4i8zkDgx/rHhxfxVX6C3uKLQRxMFy8RAEOD2z48TAyp4OD4Qr2PTg8MP2DUcWWSXQwXEHA4NrS6LojLYzBVk1yMFwxAcOETw6GKwgYW0Qo/cB5DJH5z5YAuIz4MA0GPwLGllRxMFxBwCTwFRyzIrHHLjCmMFxJfgy2wKqD4QoueTg4MPjxydt6aw4GLOwjbLk1BwMWdhLJdgcug5+9hK0+dTBgYTdhCoptMUXYT5idwbtMf7CjMDg4GLKwpzA4ugR+dhUGZwdDlsweFdWiOBiyZKykpoA5GLJkrLdoIuoS9GR8tiiOYq5VydAHHB0MWnDZG5wcDFoy9PG+OBi0ZOgDrjizRTL0AauL0FOgj92Kd/BrKdAHHBzsWgr0AScXoadAH3B2MGsp0AdcXaQ/Su8/5hBOSinQ511ANYaeAn3A4uDcUqAPOLpg9EAfMG5zXK3QB1wcbFswyRjMbggY+niPJg7M+KpjzsDEpfbxJdZFoKfm+T1rELjMMBoJMPtFyqN1ATP/gHHuEkj04eszYsn+xQL1FvfPL5ebGBZucW7RYzm7xfeL8S8aktPA1q84DfCW5jSo0+id4uKNwWmMTmNyGrPTWJxGeFRzGtVp8k6TOE3BaYpOU3Ka2DwSIwJva06TOs3eaRanOTjN0WlOTvGB5+I0wxOb06xOi3daxGkJTkt0Wgg0yU5LcVrgpc1pUafVO63itAanNTqtyWllK1qc1uq0wofVafNOmzhtwWmLTlty2tikFqfwUuPf6lS9U3zQGpxqdEqAi2anyha2OlX4vvnVscPAib1ZQjweCLOMenwxsGCP2uBhvp6AHQ/b9fBaD4P16NUe1urhp952xXBSD/v0aMkexunREzws0qMJe5ijhyN61AMPL/QwQA/X8yizHn7nYXIe3dXD3rwFLsHIPBqBh4V5+JZH/fRwLM+U82ibnsnm4UqeaeZ78A4lmFUew4yH5/hMCbPKoVR6ZLzH/OLRKT0S3qNNekS7x9rikekeS4tnR+KR6B7jikeUe1RJjwz3bBA8wttjQvHokB7R7dHUPTKbCHTxBBB4NEiPSu0R1R4F0iOjPaqjR/f0KIweLdCbaDelAbHsUZY8AtmjN3pEsVdKoDZ6JLFHYfSMOVYLlAb+EOTAmPeoJsYcQwTynj/s+BlzUxSxMohpiD0egjE3/VAY8x4YwJhbvBSmAaQqfyyQAhyMOft+xCN/wMGYm0LYI7IYc1MHhTFns47s4g8lGHM24cgjDBngYMxN+WNDjWDhDyUYcwvYMj2P7TFigT+UYMzNGCiMucXHWcSLMOYW4MY2FkbNH0ow5hZhZ7vQHhnGmJv6ZhFyprfZptAUth5Cwpibuma7M9PTbN8ltf31p59+cr9VNCcx5G+N5uwnMIgyf8uRmTlik/DtqwsybkOce8qdXQVHRz8A9453HJ/suDn0IpygOD318tny6rvV1d4pmp6wV+WcsDsY8vXqzebuxXPirwlABewv/eDpA3t/vn5+scAq0uG9+nl9/5Io92IRxss369OY8bubzV3SCdx+tv5+fb2+vLhejFnAyJu9Cj9dPlltT/aAz+COIYHB4C/Ozq5XduIGPjsn7siORvf66Xefri6ec2DID94C/wkz3xa1tpBwXGxzvg0032XZoicw+dF/mhbacP4dLfzX/zQt3A3QLxzDD5dXe+elPlxebSeFhWBzBI1lev7gq74mPrpa/tCPY3T4i1ebm6MfHZhPf3RgPgDyxavNRz3evh+aIwieRWTL6ItXm5kl0IgvXm3u2xmzOev99Xwu4CSCngyW+Gy94aTcFv768vLcAuhJ6EdX7l1ebC5fX13PpxbubmZyjjjm3c2GRWxM6md4QfiFzIC1QvPn00kzxImIwc7nfXzx7OOrq8v5iBor20DLDqr7ry+ezmyBl4B7XAxwHkLechRmzkz7AefMrHvAvfH+dPV8dfFs/7wN1PXUPQ5LRTeJW9zbw4VUsWULu4z03zwV3eLsE85CrK6PePic+tWr5VPOAhju3YG/vTbsTvvNadC4y3dIzS7rNnmb9Qi15Ttu917izWGjT9bXTMh9ekiivpmc4mn2Nt8Wce+dbdZt6pzxiBpyfba+WL98/fJfV1eXN0c9eHFwBtNYfD/18uBqdba6+p+f3uTu6Xsd1xP2mwml+6k37eypH63OPlmM2TNou5RvF2M9THm4wA+wl+XRnPBguT/5HiwP5hbId0k3mC3p9Ozpg+Wzg6bTdw+Wz44Os/bEW46zPlg+Y7I/vOmeOeXRQQqycT6rBML10+/mk0oPlq/6odGHM9PYJTxa4J5cnH319Gq1uri/fGrcB0pga3vdD8g62Ju2JO2Px7bU3vohC+DN6iFlO4H6tOp5rl7CRRdlSHaQhkTjGcqpRCBjMWxOO/jt9rwT7z4BEGvH5mr96qPV0/XL5fn17lCRseRZ2wk7lWGvdZbhqHmWtt8+Jokl7jVwB+9EVVdSUED6IeddqR0LtDK9ObrLTXu29dOcznwu/3Rxsbr6kuaRk6Vm1V4vxr9wkOjOHf5IuCNs+u+kO20LlzvpjhSg4MmQ5iw3kCXtpwd/mH33rvn5XfN3kr+T7iTLabh/oz9/hbmsls9WV0hrOzNl3baD7q8397eTJs+Txk5iMYq7NzaXrNeeLs+tMKP/vy7XFyRuFYF7y1f74Nfrlzt1sram0pJxjD+9XD5fUdGOwd9bXjw7X337Yn393erqS849dsbf0z+8fDOn9dHrqUbJ3vHNP68vz9cX29T54GLPem999fT8mNvPrzh3CtF7AvAhKvfHb1493Fd7tomP9hMf3ZZzm3iQk4yfLd98tH5u1wYwCb+42ry4vLd8ubpaztznN9yuGU/Z3mVwrOKYjHvbdg2+QQftMSzALe+Z3+5135FgPN2Ecayedckxev7fjs0Dh8Vt55H3WQfoOjyLil+9y+5/8/k9uqJf3PCz3faJfPS2Tst2EHbe5Z5x4cDjzeOcV/HsSXm2ehZZ930hLsbF/fXmzubyzrPlZnnnh/XmxZ1Xq6sPVpwLvbN8+nT1arO8eGrDszcMXb/4Bwfil++GP3+6OkceYNlZnL1lNztvZd+x939j26i9bjAV+z/1/petxZM1W3dsYIuzh7uD9w935+4f9gycBT+7v766njcyny63T9wrErZC7OXqo/X1q/Pl3sF1uOqOL9N/tuO/OUj+2eWzT5dPZvhnduTvN0o/Ho+SnVe+vvN//sed/3bH3/nvEPD7GbMB1cXk206NmI0dP2O12Bu1W8bM0v7pY/YWG8P7jdnfjsfs9zVI/0WG6Mgssmcz8bb4P77YXK1RMK3B169f/jAvah63C5zn7cUZ9jy/+Gz5ho3frfJw3hRy70i3lXx+efVyu+OA7cyGjX7hxNlXHZ0xnIOribARcFPRqT3y/e+w6LdG3OxKtrdIGLYP1xdsKj6+uvqCWy8gDfiL71dXZ+eXP2B4wb55dQV7O/S9FZxy5nszScMNDeaBO76u5Ne3cL9isbxb+H95eck+nMy3yTlzLO0u9Xjx+CX2so+WmyUFjEEtxgUVXD6/Wr68c3l2Z5fj8avzy83jx1b1vsCPdNY/IPDp/gOB32s8MICj6x0bwD9/dbm+4KqrWab93FD57ib1A+YwP+CG8gMuGT/gnvADDmY/BHuO9hztOdlzsudsz9meiz0Xe672XO252XOzZxyyfsCdKYZXDK8YXjG8YnjF8IrhFcMrhlcMrxheMbxieMXwiuEVwyuGVwyvGF4uz8kuGN5geIPhDYY3GN5geIPhDYY3GN5geIPhDYY3GN5geIPhDYY3GN5geIPhDYY3Gt5oeKPhjYY3Gt5oeKPhjYY3Gt5oeKPhjYY3Gt5oeKPhjYY3Gt5oeKPhjYY3Gd5keJPhTYY3Gd5keJPhTYY3Gd5keJPhTYY3Gd5keJPhTYY3Gd5keJPhTeC1W3Mevcf8qym70pIrQjhJcDk0l4O4VKpLqbio1cWWXUyEQ4iLnrCF7ALxuBZGbbG0Fu5jXjZ8c7h3cRV23yJ+OdfqYBdvbf85TUcJ5ThHjdvrxOxfJNa+6FGpelJNO8qRykE1Fs5PfHDlrEDeXVlm/4ILhwiaC+2guLh4nBC4rG3vn4v+uMhhjnaa46QO0aM6LGj2CJEcNo3I4+McuEudHOJXd9gE8zofFTxE7+NpwlE/Z1cPKPbBHQ1ocYfk0tuHRcQdtfCWIvGQ0HZaR3pnpYcTi347rFSd8N9x4knF785x0sAT1EdD0+zqkuMWQM5tGe3akuPMt2W0S30O29OKMYlfQ7XZSefdFV0nitdOS9tqYgf61cdc3fYe3GonLU2O/T7+pr2lr4oofytx/7ysnZ31f0wdfxvu28h5r4K/Zq7fmIg/qoe//cdMgH9+378LY/oFPXFLXbckHSyF33f1vx2p/Y7Ghy/Wz1/8wcT3JMBt82Xv9f7zbVnfixf/3QXfr/r3y/UbE/FH9X8w8Zt/v28u+665+g9Sf7Ae/8G6DkntTJw7k9+Dh4c6hNZylVpzC9xKVgZJ1XvVUFvl7EEaWgy1hJJbyFJIiCFJqzEVJQA5YJnx2oIWL7VyRHxo0gi85hptjmgMxYeUk685cBF1kIFQaB85BYTVS3RQKammVqVYzLQORYL6oEFjICK8Daml4GvMUbmSnIRcJdbYSkgEFtchRq2iRVK2sO4y5Bp9yqWGkgmgz0NsUVMuGkomVhq7DSHMPkojWj0NObUaS9HcGnHtaUjEp0st2a68jEMpIZcsooXTCBIH4vBDbMTkE+wdhiqpVi0hFgveDrQ15xJTCdXiuofiU05VS1KLVZchcC4gRI0aLEB+yKVKKy1V+k/8UHOqrWbVkDhHMcQStOYQW/ZEdOtQgxYRjZDlGkap1qJvXiXH6NqQQ4qthZYjPebqULJGCYGBJy6+DTGVXEqsoonzBG0oyauPSWpULH5tCK2klFKsQasEV4aaa25eNVUlIL4OqXj13msOSvh2GYI2eF6thTHLg5RYK92Vc25YBVokqDsGQvcrpgaRqjFoCpl74V0eUio5hJRzisWTELi22zetGlOhDolc4004fmOqMIi5+cL18CHGcFSHiMtDjjWXHENuPtdjLJxQOq6DSe+DaKyp5ITpZ6AFJcToU6yZIr4FabGpxkTr46C5Ns/5i+Rji2Zb1Ni815Zr9iSEknPVRJcRo/+OBA7pD7W02IrPhflEHUlyKFk1Rc/RkTBoiqlwnCNk5te7EhKm2MaqT7n4Urgo4BCL2XBVo8SqgUZhWObefQqw0sF6UAUG4nhcRCXZdf0lhMCJi+M65DAHjT2qVIYSQ/WptMJhjJM6gvNDCzWkEoLXWjJXz5ykHCEOp4h/cY54StoRXq43Pa73NuKOip02wCq6Jdc7K+8FMSEd9khss8h4T70/tCH6or5Jq7VyB0MdpBRJNYSiieuh89BgAVFa8JUJloeYUs2sDGHmB8x/gSVfS6vcjBAHQqVii9kH5nUwz0CR1jRELtgIg8SUfc1VlINRmFqTFG/M2SatH0psSVoIJVdOAumQOQ7jNSdfOIKlJkSyL9qQWAiRWFpuPtrq5rz+UDziLanXxNGZMiQkUw6sJ04WZVofWHMxVM4C5SG3KDFKTHZqBpYFm601CWZ3SUNRFZVcJDcO16RBko+StRU7yhTh7xJTilIyh4XiIFlaDa35xIE4CUNBFiTfUrYzQGEIme9TMAiFw00yIPeChlgSZ+dEhuYl+pxDrM3OY9GSqKFFX5OduRoqYstoCRylg8OXrHDzVDlTR8t8SjEW1Vira0MNKZfUGMzqkUTJl1RqbjFK4UjdQN9WEa1I5oKsSqXWVpKXUH10FZ6fiI1LnGQqJpxqq5qSlBCaujrEWkPxNbeqhYuMh4ACEVsqOYeasEj7FJoGX3LgiFkbpDUJUmMNCaWkDDk2n0qsoZTKGd8h5pRS8RK8R5DUIUhLLebmkcXZ5aHkLE20ijSOlh7Voa4MJVRqiK22LnsOsOhpHXlIMaOTFEkoRbh+tOYaIcyrpwhTSzQ2bS02lwffoi8502M4A9OgudXYUgw1cGwOpcTaHaLwmZJ3JSCJUpYk6DAFnGnIWZuozzWqhIgzCwUmSayFz7i8MwHBU6sPzdfgk3mUDpFkJGSLJUlOWRsaWhgKUxHVKGVcnodV+HpaJA5SQs1R4C0cxLuljoMct1QaBkRVyy2G1qD0qA4qDaoI0awcVCu3JRyQTuuP0J5U+q4ciMgTwg7oQDAJel9KKF6ctD3CQh0HRW4l/SDHrZUe5Ci/HO0tlZ6QfpJwRPrfVcdxpUiyxz8SNGQhYovxluCtx7uPwzzuH4SpNR5INLd4fHF59fLPy3O+d+ENfLm+WJ5/uL6Ygz77SbjHF1+tn79cWqTG49Uu+IUtrls8vlr+sAuI+UDc4vHl1fo51XxrX/mwGJT9w9f7ASBzVME/zavL2WouwK/mna2eO8iKOi5bK/teVmen353d2sRz7Clca2CeTc6a2wXU3V8prjkJTu2/6porLrjqiksuuOSqiy6YUwyo/2J3hnXPHSoTt/XZX7sS9zc44Hn+3uEv915fff/WoOEcctwPgbXolscEKv1l89f9GJgHV5d8PG19ebELguH1fuCLRa68PfCFQPSfDTk+CXwhIOsg7GUvYXfucy/shaDCnw17KT3s5QMWAvteezDd0mMOxVUhPKCRm4EUXTeSgrEj8WA//hR+PFR+PDR+POCH49A1ES4Eu/CAmjbALTj3PKB5E6MwJH482I8/MDLOkQ/sEyxWho0Lp74H5Zedea5hoggFdMqBlL57YQdDCtsjLhIb7McfmGzhwTZKPOC95oa9ASGHDMI/jCgSfjzY3o4HXkUeEj8e7McfBHXhAUFYeWATSdDGoLabdIh67t8Y0COEh8CPB2R45MHUAR7sxx+keeGh8uOh8eMBVUCzs/+omT/ev0fYSXE3EWp/qaUOMht4uP/hBBY/ZB8y6m5NaGe5DT6zUQ5eS0queFRzbEYSI8fsc0kDtg01bcwCWsLQYpamnJ0nsqWlISSzauQcg0upUUJVBdVbXUJV11prq7U0i39JQ43VY7DBVuViLtiMolTfEhdlREI+IjvsqDVUQkDyUGqN6vneXQgutDhw3j2ZNQzdoQQUyRAztTKHUhhawcqiLUX2OiENGlG1c4qV6z98G2LM3pfYMvdfiDJraim55Zbsgsk2VCw3pXrixKRkLDiBLX5OXBKRTWRVDVFb5D6ARBADNwWUIpkbBqIMmSsUovKNNvYQ7AhsW9K89nsEBrq7hKQx2l0VaUitZFTi3Ljvo/CVtyg1Jy9Vm2s6FPY5kaZgbGhx8JKzFrMCpexqGVprUqTUGEpxFeNWsFHE5iWu5EGxyHmvSWLjrs3Bt5YkhJw8d1Vk7HZKV0XlRpAchuTRfc3YlIpjoKP9V1LgTqHEJA7Kdoz9X3EJQ0HSELNP2kog2Ie9ZMAkllJQFzMdnLOX6jN3ZURiytQ6w1eLlGI7WWuOWXJqGD7bEGuJRULJpXKXYBlCE5R2Bp7PEA1MaL67VxL30oSAGSn6EFWwWrBntW0K2dVMksoE1BK0FB/tgrgBo2qprdUUMP6VISuvJJZqNzDkIZZcNVafsZ6ynWQj0G0+2HDZPjaR1jI6M1dYhCH4FmrJSbJdiiFD9FUq5pVml5YMqYbc0H9DNINiUc1aMUq25LEPaotsQBWzYmArGIPPNfmWC8Y0V4fqSwoVc14JjY0euwku3CjqiWcsQ+FyDg3eh9A89kEJocWIbaAkTWyfKkZH30rf+GLtw8SGuc++1ZeGVswgGJO2YJuh7HNuGBswesHQpEYfWvXFZ8nw0layN6OsZ/kRAKg15pxrDSmSIUrOEpRdKDePMBV8Zadc2SUhClpKEnNJtWG8ND08Bd+k+hgLa3pI8JqcayrK/SmIES2K+bQK3REGr0H5qGTJgpVEBo0thYIFsqI5DVVzTJpbiT6yhcfY3GIsjahnbnEZWCLaGobgysAOKZQUotRUSmb6DDFwqYcWvrzJTS9DYDgidoHQuGBkkBixAGTMqtgFBs83PWv0uQlbPwJLA8Y9qpZSTQZ7n/iOpCSPZRrxHKMWQ5aEu338QG/62mhcC2qyXHLUrDFgiBTyYD3X4FVjkeSppwgWmJollBJNQcgYODVI4M4y4qiGHArG36gZ244pEK3CI3K11qNtpCxZvdq0aWop7OLZ76aYoukXkW10qan5Ijl1jaRhRm4xa258VWqI70oxY9nvIdrqcJ9zu1Z5utc5EtWHW51f/WjZ7zW0HPX37fo1b/9B/fqPsPI/wsr/C4SViye62u7csous7I5juwwKtRXlw3yzXPKKNR1PLlcAcnEulwlz/SNuVUzaXPbOTall/vKj3QLFzWHcaYaqeRh/WY+Cck/DyOtRXLKplvv/uLi9HNdbDiN3gyvcf4iNI6mLh6HBfN/ngKzq4mFANd8AO4jbxR5yUKS5eBjZy9d1DnIoiuRNWAQhIuEQreJiOcxxGld+XKnk4yInQeV4CVAO9zuNG7rYFO7/q6dR44cUmydgv8RtYeRHI4Hhf79N5kQ8SLiFjsM2eXbVB2jRKQ/quCU23Z/WctTdGOWPasFAcITJnyI/nD4EbJ4Gih/Vg+p4gIsY8ZMkc9vt/t7++og+AtDNgPa7UVj+CA8/CRo8DN05eb1f4Las+xPnPyo6+HdBxN/dOX9Q/zMhib/XznnXaP+q8XR9Ue6vtN939b8dqT1M5I/w8OOjQbdNx5NMPeW2rPtT6w8mvvt3uu7e2jm/iy78u4f2/1Pq39Vfv28u+xtTfzAnftWe6Ez8PcPDiW5uRVOKGGwLYVVHAbtHccNs5ocQCZn2WmK1S7TxLUggbloidu+AvT3g6Qoef45wmUCt3mMcNzs3zlhCyUtsWrEshCETYkiccGyYEsJQNUorNcTGNfJCtJFKzRKrhdCFgQgqLUXx7FlM3WGGOBBd56uPKZgrJg0SCUnHyN1oSRwiLrlaMGcbTvWqPuWADdmCv7NkotzwkREDGgZfVFSJp7EEwe0UG32Rsl2fPVQcnHgqUuNCcyGyLxNSnEoxrxyRfbkWranZSXw/KP4jybGkMsftHYSDi8ffU4InQqrf2D4Q0ShRM4F74pRhaznFnHFosfW2UDqNOJK4gfwkQrwNTULFlRdq5A7wNnBhvA8FnwChkXiqCNKn+3y2sHMpqdWQQyLu0g6O7+cgsg+Pi1RfWyB0mLPmqeLhwPnE3QXHIeNl0Fp8wCMStCgx5A3/bUteY8aNUgdptcVcfcJ1S6Qf3gsvlVvV8YeloRJnSff6yvcabon/JjY/EmcmxFaz/y/FS9IgOIosKm8voeBtLzllXHuEFRphB8HtuNhPQsb3A8LxLg9tPyDcXFX7CcQK4C6y8wktJCHhJGR8Pw49ETFwlIPFpEW5yB3vza0B4Uex2kfrHufVQSQ2dRwkWDTEYWD6UQw1dZzQcZBgMRmHIeK3BJEfxZ3742qp5SAPFqAjPEcB0Jz8OAwsN9vJYTQ6vrUD6tTiuU/i0Q8ScLcdEHNbxPdplp2hx0KeTmog+TiE+zTNvrp6lC/Hmf2/Z6h3GgIMVHHjer7mcRR4S4jwUQAwp1mqZgI5CtyCBOJOtZZUOKFhfDXnxj/u45EI04zeE6ZgHzSIg2jQpFIJZ+58NpesnmhkOzkQOZNELIRw2oLoB65VUTu0lAM2zDhkeG4VDQGn5mmOZEeOclMCJDCFcjolc1TIe8l8qCENyScfUmmJL5gIEa4we81wMT6jMuQcEp8mqMGCsIkN8VRSfeE4C3HbPpVsR6kSQUl2u01NPnLeyWzhQ1GfJXh85DjvLSzSI+m0Ii1l0JZwSzcYDlb448BuQeLWHFL1+KjJQQhx5fMVuOQREK2FGj2Hn1oP7ObMlHD6JfLRkZNIbx1a0Rx9az4nvLBK7xATHwjSbgQEVA2EpoZSuXJOhxAKXyvhyxc1yHEObiupTX2IfNPDDnW1gUjrkIovqcZEyMBR5HdDPpTIgZZS7FNMFl7ifY4+tkwEyhByjQjMoB7VoHLoiDDwzEGhwOGeVnLgtFomDIjLSU7CuKsgGJP3nnoITCgphuQL0iAQurCfgBSqOK5L8onvalDiIETdjhQdBn4fxXUfhXWbRDmJ82600tecS2vEuZwGfh8EkxPRdRQaHofYQmJWByEk4DSw+yRM+2jRWzDaUaT3Pmxi6yi+/CRM+YSMwwRE4UE8NYLuXdHjR7G/hNrtt5UQjJOg7eMMR4HPhxHZuKuOA7BvofOAcIILT2K2T8K89xMseOSAUITpCWFHOU6wHCUcw8R5HNH1ixOOO/SWtp7UedSSW8g4ynFr9PlhbyDDDoMfTm/oPAl8CIf8nXjsfyTImw3arxPj/c9zsVa+Ol2jq43PVqMp27etDx2gLgdH+FK02O7g+Iow33Lmoi2CvHt4N1/tNOchn+Jx2anDleuqBXkHiydNFuztHVG2RPDyhfjkEi6xOew7HOg5hHmbvnPjt/rVw1HuPvzTV1/d/ezjt0d8/8vV6mwxSv7VUT85//yj971l+cHye7tRHjJvu2qRD/3shYv/I8EsRIpzOfUHXKBH9Kna+Vupqw8IYeG+6pbKUDgEJpr5+FW/vToNZf6nQtjtfJm1EtuqnMNGcgl1PBRazi37duP1FrJTsr6lqEH5KCDVzhlbvxp7hlS7NSSxwahcCbl/PTYN+NMFH0ehiq9eLJ9d/rB3b/e9yytuwV8+W7++7vdx7i6s2o3IfGf1hhtzLVD/rd884SZLXvaPHwVoOb75f/fxo5tPE82X0H96efF8xZXHsfZo/OfrC7qlf65ofcHFqUcfYvm5W0mPsy43qzdvmzAQPrfSGmgt/QW3cL9vj2w7aP/zCLsO4eX8iQm+YPBL5nC/B/TletPvfO1HFjgF8YUdqdndfJxitJ1sDLH5oD/91S0uuYf1Lxevz885OsIV1F89vVq/4viN9cUn6+cvztfPX2zuXV5crJ5ubr6KcX/9ZvWsV322PL9e/fT/ALkH/s0=').then(json => {\n",
       "   const obj = Core.parse(json);\n",
       "   Core.draw('root_plot_1779222679379', 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_1779222679379();\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
}
