{
 "cells": [
  {
   "cell_type": "markdown",
   "id": "ec9293c3",
   "metadata": {},
   "source": [
    "# multifit\n",
    " Fitting multiple functions to different ranges of a 1-D histogram\n",
    "     Example showing how to fit in a sub-range of an histogram\n",
    " A histogram is created and filled with the bin contents and errors\n",
    " defined in the table below.\n",
    " Three Gaussians are fitted in sub-ranges of this histogram.\n",
    " A new function (a sum of 3 Gaussians) is fitted on another subrange\n",
    " Note that when fitting simple functions, such as Gaussians, the initial\n",
    " values of parameters are automatically computed by ROOT.\n",
    " In the more complicated case of the sum of 3 Gaussians, the initial values\n",
    " of parameters must be given. In this particular case, the initial values\n",
    " are taken from the result of the individual fits.\n",
    "\n",
    "\n",
    "\n",
    "\n",
    "**Author:**  Jonas Rembser, Rene Brun (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:25 PM.</small></i>"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 1,
   "id": "d363fda9",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:25:13.883983Z",
     "iopub.status.busy": "2026-05-19T20:25:13.883863Z",
     "iopub.status.idle": "2026-05-19T20:25:14.952117Z",
     "shell.execute_reply": "2026-05-19T20:25:14.951488Z"
    }
   },
   "outputs": [],
   "source": [
    "import ROOT\n",
    "\n",
    "import numpy as np\n",
    "\n",
    "n_x = 49"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "7764e091",
   "metadata": {},
   "source": [
    "fmt: off"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 2,
   "id": "008d347c",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:25:14.954453Z",
     "iopub.status.busy": "2026-05-19T20:25:14.954253Z",
     "iopub.status.idle": "2026-05-19T20:25:15.070472Z",
     "shell.execute_reply": "2026-05-19T20:25:15.068712Z"
    }
   },
   "outputs": [],
   "source": [
    "x = np.array( [ 1.913521, 1.953769, 2.347435, 2.883654, 3.493567, 4.047560,\n",
    "    4.337210, 4.364347, 4.563004, 5.054247, 5.194183, 5.380521, 5.303213,\n",
    "    5.384578, 5.563983, 5.728500, 5.685752, 5.080029, 4.251809, 3.372246,\n",
    "    2.207432, 1.227541, 0.8597788, 0.8220503, 0.8046592, 0.7684097, 0.7469761,\n",
    "    0.8019787, 0.8362375, 0.8744895, 0.9143721, 0.9462768, 0.9285364,\n",
    "    0.8954604, 0.8410891, 0.7853871, 0.7100883, 0.6938808, 0.7363682,\n",
    "    0.7032954, 0.6029015, 0.5600163, 0.7477068, 1.188785, 1.938228, 2.602717,\n",
    "    3.472962, 4.465014, 5.177035, ], dtype=np.float32,)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "d13ea803",
   "metadata": {},
   "source": [
    "fmt: on"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "33f6e1bd",
   "metadata": {},
   "source": [
    "The histogram are filled with bins defined in the array x."
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 3,
   "id": "521727d2",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:25:15.072674Z",
     "iopub.status.busy": "2026-05-19T20:25:15.072519Z",
     "iopub.status.idle": "2026-05-19T20:25:15.203481Z",
     "shell.execute_reply": "2026-05-19T20:25:15.202698Z"
    }
   },
   "outputs": [],
   "source": [
    "h = ROOT.TH1F(\"h\", \"Example of several fits in subranges\", n_x, 85, 134)\n",
    "h.SetMaximum(7)\n",
    "\n",
    "for i, x_i in enumerate(x):\n",
    "    h.SetBinContent(i + 1, x[i])"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "4d344418",
   "metadata": {},
   "source": [
    "Define the parameter array for the total function."
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 4,
   "id": "dea638bf",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:25:15.205434Z",
     "iopub.status.busy": "2026-05-19T20:25:15.205283Z",
     "iopub.status.idle": "2026-05-19T20:25:15.308880Z",
     "shell.execute_reply": "2026-05-19T20:25:15.308115Z"
    }
   },
   "outputs": [],
   "source": [
    "par = np.zeros(9)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "e3ce7fb9",
   "metadata": {},
   "source": [
    "Three TF1 objects are created, one for each subrange."
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 5,
   "id": "760568f7",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:25:15.311043Z",
     "iopub.status.busy": "2026-05-19T20:25:15.310910Z",
     "iopub.status.idle": "2026-05-19T20:25:15.438011Z",
     "shell.execute_reply": "2026-05-19T20:25:15.437205Z"
    }
   },
   "outputs": [],
   "source": [
    "g1 = ROOT.TF1(\"g1\", \"gaus\", 85, 95)\n",
    "g2 = ROOT.TF1(\"g2\", \"gaus\", 98, 108)\n",
    "g3 = ROOT.TF1(\"g3\", \"gaus\", 110, 121)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "a5d1a539",
   "metadata": {},
   "source": [
    "The total is the sum of the three, each has three parameters."
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 6,
   "id": "2d432f90",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:25:15.439986Z",
     "iopub.status.busy": "2026-05-19T20:25:15.439850Z",
     "iopub.status.idle": "2026-05-19T20:25:15.554735Z",
     "shell.execute_reply": "2026-05-19T20:25:15.553996Z"
    }
   },
   "outputs": [],
   "source": [
    "total = ROOT.TF1(\"total\", \"gaus(0)+gaus(3)+gaus(6)\", 85, 125)\n",
    "total.SetLineColor(2)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "bef091f6",
   "metadata": {},
   "source": [
    "The canvas that the histograms and fit functions are drawn on."
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 7,
   "id": "7770af6b",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:25:15.556863Z",
     "iopub.status.busy": "2026-05-19T20:25:15.556723Z",
     "iopub.status.idle": "2026-05-19T20:25:15.694727Z",
     "shell.execute_reply": "2026-05-19T20:25:15.694009Z"
    }
   },
   "outputs": [],
   "source": [
    "c = ROOT.TCanvas(\"multifit\", \"multifit\", 800, 400)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "140c48a1",
   "metadata": {},
   "source": [
    "Fit each function and add it to the list of functions. By default, TH1::Fit()\n",
    "fits the function on the defined histogram range. You can specify the \"R\"\n",
    "option in the second parameter of TH1::Fit() to restrict the fit to the range\n",
    "specified in the TF1 constructor. Alternatively, you can also specify the\n",
    "range in the call to TH1::Fit(), which we demonstrate here with the 3rd\n",
    "Gaussian. The \"+\" option needs to be added to the later fits to not replace\n",
    "existing fitted functions in the histogram."
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 8,
   "id": "ec0c9c92",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:25:15.696875Z",
     "iopub.status.busy": "2026-05-19T20:25:15.696735Z",
     "iopub.status.idle": "2026-05-19T20:25:15.841767Z",
     "shell.execute_reply": "2026-05-19T20:25:15.841050Z"
    }
   },
   "outputs": [
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "****************************************\n",
      "Minimizer is Minuit2 / Migrad\n",
      "Chi2                      =    0.0848003\n",
      "NDf                       =            7\n",
      "Edm                       =  8.86911e-08\n",
      "NCalls                    =          106\n",
      "Constant                  =      4.96664   +/-   2.83221     \n",
      "Mean                      =      95.4663   +/-   12.3905     \n",
      "Sigma                     =      6.82779   +/-   7.49131      \t (limited)\n",
      "****************************************\n",
      "Minimizer is Minuit2 / Migrad\n",
      "Chi2                      =    0.0771026\n",
      "NDf                       =            7\n",
      "Edm                       =  1.00182e-07\n",
      "NCalls                    =           73\n",
      "Constant                  =      5.96312   +/-   1.14355     \n",
      "Mean                      =      100.467   +/-   1.53372     \n",
      "Sigma                     =      3.54806   +/-   1.16899      \t (limited)\n",
      "****************************************\n",
      "Minimizer is Minuit2 / Migrad\n",
      "Chi2                      =   0.00877492\n",
      "NDf                       =            8\n",
      "Edm                       =  4.98832e-06\n",
      "NCalls                    =           87\n",
      "Constant                  =     0.912053   +/-   0.435309    \n",
      "Mean                      =      116.304   +/-   8.32344     \n",
      "Sigma                     =      8.38103   +/-   18.5139      \t (limited)\n"
     ]
    }
   ],
   "source": [
    "h.Fit(g1, \"R\")\n",
    "h.Fit(g2, \"R+\")\n",
    "h.Fit(g3, \"+\", \"\", 110, 121);"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "8c3eaed5",
   "metadata": {},
   "source": [
    "Get the parameters from the fit."
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 9,
   "id": "0929ed36",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:25:15.844759Z",
     "iopub.status.busy": "2026-05-19T20:25:15.844619Z",
     "iopub.status.idle": "2026-05-19T20:25:15.953548Z",
     "shell.execute_reply": "2026-05-19T20:25:15.952796Z"
    }
   },
   "outputs": [
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "[  4.96663958  95.46632975   6.8277931    5.9631179  100.46745499\n",
      "   3.54806038   0.91205321 116.30403822   8.3810307 ]\n"
     ]
    }
   ],
   "source": [
    "g1.GetParameters(par[:3])\n",
    "g2.GetParameters(par[3:6])\n",
    "g3.GetParameters(par[6:])\n",
    "\n",
    "print(par)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "fb4ac6bc",
   "metadata": {},
   "source": [
    "Use the parameters on the sum."
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 10,
   "id": "5088a0e1",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:25:15.955265Z",
     "iopub.status.busy": "2026-05-19T20:25:15.955134Z",
     "iopub.status.idle": "2026-05-19T20:25:16.075849Z",
     "shell.execute_reply": "2026-05-19T20:25:16.075124Z"
    }
   },
   "outputs": [
    {
     "data": {
      "text/plain": [
       "<cppyy.gbl.TFitResultPtr object at 0x556536e2a390>"
      ]
     },
     "execution_count": 10,
     "metadata": {},
     "output_type": "execute_result"
    },
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "****************************************\n",
      "Minimizer is Minuit2 / Migrad\n",
      "Chi2                      =      0.31282\n",
      "NDf                       =           31\n",
      "Edm                       =  3.25007e-06\n",
      "NCalls                    =          495\n",
      "p0                        =      4.91052   +/-   1.41324     \n",
      "p1                        =      94.4492   +/-   3.71244     \n",
      "p2                        =       5.9461   +/-   2.41662     \n",
      "p3                        =      3.22456   +/-   3.11384     \n",
      "p4                        =      101.662   +/-   1.67862     \n",
      "p5                        =      2.48631   +/-   1.91151     \n",
      "p6                        =     0.911626   +/-   0.368736    \n",
      "p7                        =      117.581   +/-   5.06092     \n",
      "p8                        =      7.59194   +/-   8.78217     \n"
     ]
    },
    {
     "name": "stderr",
     "output_type": "stream",
     "text": [
      "Error in <TFitResultPtr>: TFitResult is empty - use the fit option S\n",
      "Error in <TFitResultPtr>: TFitResult is empty - use the fit option S\n",
      "Error in <TFitResultPtr>: TFitResult is empty - use the fit option S\n",
      "Error in <TFitResultPtr>: TFitResult is empty - use the fit option S\n",
      "Error in <TFitResultPtr>: TFitResult is empty - use the fit option S\n",
      "Error in <TFitResultPtr>: TFitResult is empty - use the fit option S\n",
      "Error in <TFitResultPtr>: TFitResult is empty - use the fit option S\n",
      "Error in <TFitResultPtr>: TFitResult is empty - use the fit option S\n",
      "Error in <TFitResultPtr>: TFitResult is empty - use the fit option S\n",
      "Error in <TFitResultPtr>: TFitResult is empty - use the fit option S\n",
      "Error in <TFitResultPtr>: TFitResult is empty - use the fit option S\n",
      "Error in <TFitResultPtr>: TFitResult is empty - use the fit option S\n",
      "Error in <TFitResultPtr>: TFitResult is empty - use the fit option S\n",
      "Error in <TFitResultPtr>: TFitResult is empty - use the fit option S\n",
      "Error in <TFitResultPtr>: TFitResult is empty - use the fit option S\n",
      "Error in <TFitResultPtr>: TFitResult is empty - use the fit option S\n",
      "Error in <TFitResultPtr>: TFitResult is empty - use the fit option S\n",
      "Error in <TFitResultPtr>: TFitResult is empty - use the fit option S\n",
      "Error in <TFitResultPtr>: TFitResult is empty - use the fit option S\n",
      "Error in <TFitResultPtr>: TFitResult is empty - use the fit option S\n",
      "Error in <TFitResultPtr>: TFitResult is empty - use the fit option S\n",
      "Error in <TFitResultPtr>: TFitResult is empty - use the fit option S\n"
     ]
    }
   ],
   "source": [
    "total.SetParameters(par)\n",
    "\n",
    "h.Draw()\n",
    "h.Fit(total, \"R+\")"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "57ca1f02",
   "metadata": {},
   "source": [
    "Save the plot for later inspection."
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 11,
   "id": "dc34c729",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:25:16.077625Z",
     "iopub.status.busy": "2026-05-19T20:25:16.077478Z",
     "iopub.status.idle": "2026-05-19T20:25:16.234944Z",
     "shell.execute_reply": "2026-05-19T20:25:16.234188Z"
    }
   },
   "outputs": [
    {
     "name": "stderr",
     "output_type": "stream",
     "text": [
      "Info in <TCanvas::Print>: png file multifit.png has been created\n"
     ]
    }
   ],
   "source": [
    "c.SaveAs(\"multifit.png\")"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "547f9c08",
   "metadata": {},
   "source": [
    "Draw all canvases "
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 12,
   "id": "e1399fac",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:25:16.236676Z",
     "iopub.status.busy": "2026-05-19T20:25:16.236517Z",
     "iopub.status.idle": "2026-05-19T20:25:16.456728Z",
     "shell.execute_reply": "2026-05-19T20:25:16.456035Z"
    }
   },
   "outputs": [
    {
     "data": {
      "text/html": [
       "\n",
       "\n",
       "<div id=\"root_plot_1779222316447\" style=\"width: 800px; height: 400px; position: relative\">\n",
       "</div>\n",
       "\n",
       "</div>\n",
       "<script>\n",
       "   function process_root_plot_1779222316447() {\n",
       "      function execCode(Core) {\n",
       "         Core.settings.HandleKeys = false;\n",
       "         \n",
       "Core.unzipJSON(34611,'WkwInjIAM4cAeAHtnXuTHLeR4L8Ko2P/sHfLdYnEM6tuL4KixJVv9WBYkkWtznHRJHuoPg1nZmeaEmWHv/vFL1HV00OKWsmyvXbsOjxUJQoFJIBEvoH+w+b/Hr672l1sX+w20+bTB9uLb7Y3n++efHKxvbr56vKwGTZnn13s//3l7tfvbiYZNmfv7A83/enjJ/9v9/RA+YZqH18d9pcXC/Cv+4tnmykOm7NjS9Mfvq+vt3UQY9KY0rA5+2B/sXtweX55vZnCAn5y+O58dwt+vn92+KqDD/fn50tlkAVcK4v457uzw4fb6+f7i80kIyW/2T//6rWidy4Ph8sXd6t9enl1t+Dx2R4kdNicfXH7eL8/0vDjm8P2QC9m1LkD3e8QHz+83r7YvY43Za8N/Fjv7oCOVdfitUk+P50Zmnzn8vrZ7vqT/e+X2Tsp/PDy2a6v6+OwmWobW01mllMsOfgAwmb6lYxFqkZpsdVqKplx6mYKScagWURSLjE1BqybqY4tZU3VtFhOwnAfHy7vP7l5tH+1O/96M/0qiOZRcpQkEix6e4fLH/N6MwUdzbKFFmosppVO77QeUx615hgkaqjmS3Xb+A+83Uy/SmmUIDW3KFKSN/7Zncahr89uWzsFN1O1MmzOfnvng1hB4Le3n9wt2Ey/6gXr9BwuH3/9vSvhs/MfvoXypBZLVWONWrNUluW29S9o/c0FWlrn7fcu9/reX4tq1VBDDUFDK07255fffvTug05LX5wCjz+78hfU+uLk+fNj6fvHp/tPbu40dP/JzZ227j+5uf3s/pOb2y8/e/WCzd0gpc++82cZNTYJVVoLNTWW4bNXL7avNlOIsJjPvnMAAvr0q91hu5kiy/noq/3ydP/mavf08JvtYX/Zx/XRyxdPdtf9+dP9069f3T5+1x8/uHy+FH5w+fy27Pf97aPts0fb/QXsYdicPbi+vLn5artfGjyCjy4XVnu6b8G+w7db9sPLZ/uz/e7ZZjrbnt/shs3Zv1zvn726C353C95/cvPg8vL6pP57z/aH7RMY2uH6JQ083L/aPbsz7rXpR9f7F/vD/pvdzRs8/YP9DSJjFScLuL2+3kxf/m7YXF4dePjjsDl779Xu6c1munh5fj5szj7q8ufFy/PD/mzvUufT/QFs7pR99PLFo+357nBYBQKT99Hu1eHN0nd//cmjD+5/sZk2/7A+Dpuzdy9fPjnfvfPy7GxdvN/sDtv9BTO3jPvxzf73u89u1vdf3AX97W922/PNpHTurzscgD/fXzy7/PbTyyu23yn8xSm8sOXbCu/vEEILMXy7so8HX20WHvFgezi8Mdn3D4curxnZ43d2h293u4tFIN2BfCofXl+++PTyajOFEQp6/Gx7gN878MUKILDvdyD8cdh8/eHlN7uPr7b//vJIF1//ZseM3C08e3///KsPGMIie51Gt4enX63T+vUnX11++943u4vDJ4ft4eXNkRS/vv/ycAkxHGt+uLt4+c72usOQyv2nENvxi7Pf7LbPPr44/2794uzz/eGry5eHU7pcafX97c1CaWvJaa0vX1NJ/myKD0zlrYrP57snvvf3F8/fpv1AGQ/Otzc3y8agXle3TguuINKNTJrzsPzNYZJBBpnVS3mKk/hbmZO/05znfPxG5rLUpLwudXluU7NBQxhamo3nlgYNZQ50mIblbw5hCqZDaDqEqnPQqZbB/z+HOAXRYfmbQ5qC1mH5m0OeQo7D8jeHMoXahuVvDnVSScPyN4c2qdqw/M3BJk06LH+zilcOxvcya7gL6hQsDsHCEEqbNU6hlSE0wDhr6v2WPIQss+YptMibIeQ8K1gtSCaZtU4h5yEkHUKUWdsUYvYBtjKrTaHmQaWAyRxlClEH77rEOYYpxAoGQygyR51ClCHkNtBwjL2f1vq3TJX2fnKcY55CqIN/n/IcyxQkgYGjHOsUxAb/IOoc2xQ0+byq2BztdtZznpNMocrgw7UwpzCxxI5zzHPSSTUMPj+ATFXrSEaZU5p8ilk1TXPKk88paASbU5loxtJgZU61d8NEpTqnNnmPwdGck00h1cEHEXTO4kg4ic05TCGVwWl2zurzsABxsrY+p6mm9TlPOazPZdK6Ptcpro8Nwl6asUmGBAHnucgkQ2n9ma1joT/rxKB8W82FHRSi9TdsIqjDP88ArXag0EfonRTfSXHpxXv3TRrmQvcOaJ0r/Tsgaa4gABCazBUMHMhlruseZn0rGPgbaXMFA4CW5+oI5DyUMNd1J8c617X7EOdqyCkqydxk6gjyHKbsw+BZp+bTwDNbuHMXmRtU6dMA0Lfv0hTbpKwAW/fYCRvXp4Fv2LY+24PMtjCu2ABC50ghAXTmFQz2ZdHfhFIAUgeSAnQuBlHJbJ2NmbdcvVYNlDd/TpVn82eNg8xBZAohDCk7UwgSnEKTDknnINqZgrG8YQ4Sp1aH4FTD6+TjCLkMETBPqm2IdchpDlKmkGWobFCbg8AyitcNhabYnjrEOKgCMiVDYVbDHJy3ZvaQziHAxsIQw6BA6rMKMyo2B3grz/xpngPMldVk4mk35I5xp+cAd5UwmA1BZA6hTjA4xhPaHOCu0F2xIcDRbQrN6DjMQWUKCZ6VhlDg7/CxPm8NUCf/TOrQbA4aJ43ej8wB1toyYmFgmwd4a4iDszKlcvGtrqF1nOGuJyJAJUwW2JIWZmVJYIn9b1aBN8Vh+ZuVNQmMx/9mldz5LWx9VhalRWcBMqvUib1Zw1BpuU0araMV6qxiE1/BrPOsLAlUwMezhjDVMPSPZw3qfHn5elYWhf74eNaQ+LJ/PGvIfNk/npX1gKOHytezhtqxXTttx8F5r+b4rt0i7mBHHelZfUnAtvfLkvSxhhZn9SXxwWpMs2o6jtZ7ZlGW0YK0wsVvZ1XrOt7+cTsOmI/juofhbZLmuG5jiAm5xU5epHiH+2b296JzZD+v7x0+0Qkc7rua+jnMkY29VHewHTUIB+FtXd3IYU5yqxjJnKTrBHwOpMcvgW71FKCEOEGizEmYjc69GF/yzey838eT2MzASOMc5wTpiB7Hl4J0mPeic0JZAkUkqMNKX8vwUogOSXL8g+NBVUYTXMjxsUMrq5c5BWf1fVwBVs8nlCNo+EDmpMiZPiRdlqTrjHNCHVqGhEhLbNplSD5k6GMZUn8PJ1Efcn+P9rEMiSnQtg6pV7d1TF47CiAdeOUYHOyrNqe4Cj3/NLrU6+s9p+iC9/ihS97jd0yHL9qcIcmFRmglWzgudYdvlx6Mst0ufn+f+nj6GOcMSfZ+HOVsBZT7qs3Z6gp5Wz72PsA5mw8dSHQuwhK4VJyLdC2db4owZl+0uQgj9kWbi3RtHY0BKuw7zFEssK++wxa4+ZCXZZsLMoX3fdnm4izMKbB3CRV2CvQhFUSL77DeXIjeHEMEQSSLE+FSuc/HOqxwaj6URWuH6/Fl3x3L4JAqUGOfDnVCYNXmAkmu4+M7ONYyPqikwLL6DlvgRfAvy1+cZ/Xx+VS7lu47bIG7eOsEUtTW4dE3WnofnUN9fy41460RVWKfFDSSgnbOHvPnPh392UmDx7oywo72XBdGuEz7XKG6BcUc5grRLSAf5+N4/eXt8jO8uqg3jsJcFwWHnrxyH11/2VYe6NyhwQOXec1hbtI5gS/K3BYm6EyiyXHNebOg03lPW6hxwXduKzU6p2nQ4kJ73sdxstExO0fuXSxESGUqIkoXtkJN54RLk7eMkDcLO+rDbcjRlSGLzg052sneK7cj5H0ceRFKrhMg3QMs0qHPztwQoStfp1mYJNTWmXGDRy4gHy+a8ILSQo3LLMxtoca1p75XgUAJm3HpCGV74V0dKYhxGRvbomEyLlPEtmixT1Knm7ndMkgwji4vwKnXXbj2IqhbXO2DXvmWZ3tHbqKsn6ZTG6Ul59i04+0m59h9880tuaW00ElaGTbjTLfogF26Rac3c0SHgaUjOr3yig4N5RUbB5YZ6iJ6bpiLy7r4QPKRjLzdDEYLd27dalwmdG65Y7RM+NwyGPlbXCHZMboFwWjhj81tyEUetNIFWl+2uXU7si/a3NyQXBhn64Zkx3ZupYv3TgpzKydClv7dnLxttc/POmq3KI+fdpNyfVcXV8iKUF2k34pS7YS9TkrthN05dqudsNeB1q4MHWehYtmsFnarnbDXicDQ7C+hqNYJu4u/1vqyrVVbX7UuNFrzNetSojXXgRZ0sDidtvriYnMuPNPnqPWtv/bRt/6xD9/6Sxe2yB7/bOHU6wwsjHrpf1EOXEy1hUkvQGeK62ouPPq4nguPXps5CiAn6IVHe0O2sOgFuOXQzJudsGgo2GDSiywWnW1h070XW/j0MhJb+fSy1Lbw6T6ZtugMK9Q1hmWlbWHUy8vFu+Y0BVKYG0x8/5sNZgkd97/Z1Fewy1DTNOWIWR0inzqZJ0zjOpuWKYtb4Amjv07JMNRDotE2pTZgeKYym9qUypBlCKnNFmVKsM8h5DBbDFNKgyOTZos6pTg4e6qzxTglHTC/C56ENCUZ8LgVnS3mKdqQDV/YbLFMsQ2luwgt1inWARO50kObYh4wgCs9GG6TUoZQy2xJphiHUodQbbYUpqhY89jGlnSKbvKFlmdLcVIbsIVbnS2lSZs7CQy3RsZtVTPewtnwpJWh1iEYPVQ8BNWGYPSAsTo0KDnMlmxSHVrX68z9aIObfnU2mKMNLeOTmg3W2IZWBg1xtoz3YGht0JBny2kKZTC8V2227AYhjjANs8EX82C6QGx5d08obaLf4XLpkJvNVjpU0KsGqwsEyQyGnY+rBp8BvggNNluJQAHHhYOOjFOtg7m/jcvb0sHuuLBSO8gg22ylMawgjBIQh+yAa8LBiutzwFWqoc5Wg4O4wxzE5TPgEuqgz1AIuoB4iQe8F7jADc7YhuC+BkC8C91f6GDtIP3yFmN+CHg5mG04IyBo5NmaOIinwMHQQdDgLb5rd3p0MHYQNNJsOOR4y/ABccgOQRk+YOkgaADWDoJGnA0GSWXQADQHI2jE2axjhe8p6GwYUQ3v4wL2ucLfFMJscEm806AB2FcwggZgX8EIGvjgnKBCBA1Apyj3rjroJBUS9G2zYdKnIeDFEnxpOO2AIfgGjDMKhrHCeNWHkHBH+Ps4BTzQZVCp1MeDiUd6hXOH2wrjlhlCMjzi1K8OZ9BxGKkC2xlUMu/dnRoy+AAHmWwRxQ4GBz2MwVudcKZ3X0mQEB0EF3+bHASVyNs84dIBEwcLYAERPI448PAnrmBzEDT8rU1485gVXIwq/q17BwEDHcHOjDlVBUnCH8aUuarr7tZgTBm6ri4wY0LbxQcJdwJtd+ItMIi6utvDHAYuuIzAHP7lyHTydx5J9xhgRthmCI3+IxrLAtM/RhihJY9KAONh7D7FBj4RX+wCg0/EOOx+UfyTEj0i4H7S7smtMN6Ae7YxG7HBlgncwMGDwN7pHxh8YPD03/l9EFg8/QPTP0ye/pEA9J8iwqL7MZmPlBAlHQaflBE0AX5e6T8VxFBozDn9p4qQ6rD33xBhPVrFfCRDwHWY+ciC+PPwWwGfHBCO7hYt9J8VQdph+l90Yo9xsR45TZnvEej0j15ch9AwHOg/lwkZ6bD3V6cCPrz3/toE7TjM+LNNhLjw+7qrvAgexA7TXwkTPmveJ/orHjDsMP2VODE3ROgS/ZU0EVp0mP5KnujLYforxR3dDjO+Uif32qNHeH9tIvgH7O74YpMR/QOmvypElTrM/CIWBIThd3SIYMDA8QJ6dJ2ZKXKOSEGaECO9gD4RD2yaZh4mDIKAwMtKgXqvMBa6pcC7JSxAtxR4t8Tt2KPiwbYgCAoltEgBK4uogHS9gG4RFh5xpAa0hbiAGIl84ucXBAbUSI1At4gMyO+2AJ88UUtqMP+IDQjMC0AMwZEWPIgWCKIDEvMaIIbwQDvzAiZo9cF5AZi6E26ZjwCm7g9Z5sM5vQsR5pSIBai7GFnm1Hk5gsQpjxpgiiiBLXkBmCJMOp8h1EgBCsGythR4EMg5QRsCbNbDQGxtLyDqskaCvMAjNMQd1jYIpSBTPBDRBvgp0SCnZmBiKx56YGhtgKF6QIiJaJVAafCAkAfQ6wCDJSIUDHqrAwyWmFAwhlUHGCxhoUAcB9ijRLhxQGiFu6fKYRD2MARDQsUDxnIBn+IRG0JDRC2pDwMmNkRI2GHwIzjkxNwIMgaiQwRZGA8MlwCROm23oYEfvhQ4PuFc+se7t5A6DNdDRHB4kwGGS5TIY87A9I93D45vYajMF/4U3whhgOESJVI4PLD3T4QRmtYBhhvw78HhgZkfnCq+TeIAww0IGEJbFoceSNNJfdekAYYbXMBA8GmA4QYEjG8ilH7gPGmE/jEJgMtEDNNh+kfA+JYqAww3IGAS/RUCsyHEnrvAmsJwQ5JJfYPVAYYbCMnDwYHpn5i877c2wIADdoRvtzbAgAOWRHLrYMj0jy3hu88GGHDAmnAFxU3MELAn3M8gQ2Y+sCgc1gGGHLApHI5Don+sCofzAEMOGRsV1aIMMOSQ8ZK6AjbAkEPGe4smYkMCn0zMFsUxeGg1ZPADjgMMOhCydzgNMOiQwY/3ZYBBhwx+wJVgdggZ/IBtiOBTwA9rRQb4dSjgB6wD7DoU8ANOQwSfAn7AeYBZhwJ+wHWIzEfp8wcNEaQMBfxkUFRj8CngBxwGOHco4AccB3V8wA+YsDmhVvADLgNsO+CScRhrCBj8eI8mDsz62gDNwMRD7etLrksAn5qX9+xB4LLAaCTA2It8j9YFDP0BE9wlkeidl2fkkv2Dp+VtHp5fbg9RN8Pm3LPHch4232ymL03TYIrpVwZya03bYGqDRRmMEG/UwWIcLKaBvF2LZbAIj2qDRRssyWApDJZ0sBQHS2mwhPFIjgi8rQ2WbLAsg+UwWNbBchwsp8GIgecyWIYntsGyDVZksBIGKzpYiYMVEk3yYKUMVuClbbBig1UZrIbBqg5W42A1DVYxRctgtQ5W4cM2WJPBWhis6WAtDtbSYA0jtQwGL3X+bYOZDEYM2nQwi4OR4GJ5MMOErYMZfN/j6vhh4MTinhAhAuGeUSEWAwsW1AaB+QoJOwLbFXitwGAFvVpgrQI/FbeK4aQC+xS0ZIFxCnqCwCIFTVhgjgJHFNQDgRcKDFDgeoIyK/A7gckJuqvA3sQTl2BkgkYgsDCBbwnqp8CxBJITtE2B2ASuJJCZ9OQdvoCqBMeMwHMk84V75VAqBRkvuF8EnVKQ8II2KYh2wdsiyHTB0yJYJIJEF5wrgigXVElBhgsGgiC8BReKoEMKolvQ1AWZLejgJDkHQYMUVGpBVAsKpCCjBdVR0D0FhVHQAsVFuysNiGVBWRIEsqA3CqJYjC9QGwVJLCiMwprjtUBp4B+SHFjzntXEmuOIQN7zDxY/a+6KIl6G4Bpiz4dgzV0/DKx5TwxgzT1fCtcAUpV/PJGCPlhz7H7EI//QB2vuCmHPyGLNXR0MrDnGOrKLf/iCNccIRx7hyKAP1tyVPwxqBAv/8AVr7glbrudhHiMW+IcvWHN3BgbW3PPjPOMlsOae4IYZC6PmH75gzT3Dzq3QnhnGmrv65hlyrre5UegKW08hYc1dXXPrzPU0t7tCbb/74x//OPylsjnJJH9rNmc/dkES+VuOySwZmyRrX19Qcc1s7iX3jg28dt4D8ORMx+vHOW4PuvjpjDdPuny4vf56d31ycqYXnDS5FBxPg3y6e3W4f/GctGsSUAH7Sxn9oIe/P98/v9jgFekVTtrn9cNLctqLZxhvX+3fzBC/fzjcp5zE7Wf7b/Y3+8uLm82UPXuaNycNfrB9sltP89Cfw72HRA8Of3x2drPzYzbw2aXwiLYfKCFB/4PdxXMOCckoJDn7Gqyf+lhIN3/9s8P5mmh+rLJ2T2LyF383I/Tl/BNG+G9/NyM8LtBPXMN3ttcnh6Te2V6vROEp2Bw7Y5ueP/qk74l3r7ff9sMXHf746nB70KMDy1mPDizHPT6+Orzb8+37QTmS4NlEvo0+vjosLIFBfHx1eOgHy5aqD/fLuYA3Muip4IXP9gdOx63wp5eX555AT0E/qPLg8uJw+fL6Zjm1cP+woPMax7x/OLCJnUn9AC9Q9tlPYAbsFYbfjyQxSCBORIwOvXfx7L3r68vlXBo720GvTlcPX148XdgCLwFPuBjgsoS85eDLUpnxAy6V2feAJ+v9we757uLZ6ekasOulJxyWhm4L177XA4U0sbKFY0WmZyHFYXP2Pmchdjev8fCl9JOr7VPOAnjfx1N+J2M4HvFbysDxWO8uNseqa/Fa9bWuvd7r4z4pvD1a9P7+BoI8xYci2lvQKX6ub623dtxnZ626ljIXy9cny0DRh/uL/YuXL/5td315e9SDF3cOXjqL76deHl3vznbX//LBbe1efjJxveB0mCzQaentOHvpu7uz9zdTFhbtWPL5Zqp3Sx5v+nGfY5UvloJH21Pie7S9Q1t0fiy67dmL3jxw+mj77M7QmbtH22dvHmB9tH32PWdYH22fQeyPb6dnKfniTgmycTmrRIf7p18vJ5Ueba/6SdHHC9M4FnyxITy5Ofvk6fVud/Fw+9S5D+jB1k6mH5B9cEK2FJ2ux/rVyf6hCuDt7qFkJaBOVr3O9Qu46KaMyQ/SUOg8wzjnCOQsph8qBfx8Pe8E8D5A8HEcrvdX7+6e7l9sz2+Oh4qcJS/ajh5VhpPReYXXhudlp+ODjrzwZIBH+CiqnJ+6eO4Hm49fHVmgv+3DsWNtxrO2z3A687n89cXF7vo3DI+abFRv9mYzfclBonv3+CfovYDRfy/daytc7qV7oQCpUCEtVW4hLzotV7lb/fiuyfKuyb0k99K95DW977/QP7+Duey2z3bXSGs/M+XTdoQe7g8PV6LJC9H4SSxW8fjGacln7en23D9m9f/35f6CwlUReLC9OgU/3b84qpO1NQstOcf49Yvt8x0NHRn8g+3Fs/Pd51/tb77eXf9me/F8OYrdy9+5fLWU9dXrpY7JyWHN3+4vz/cXa+lycLFXfbC/fnr+OrdfXnHKFKRPOO9jVO73Xl09PlV71sIvTgu/+L6aa+GdmlT8cPvq3f1zvyoAIvz4+vDV5YPti931duE+f0FzzXnKen/B6yqOy7i3mWvwDSbohGEBrrxneXsyfa8JxjeNMM7S+0lkTs9/zzFkPzvvZ5D9uPzxaG9nDfTYuckiLf6zZu398PCH5mwxcb9yxRE+tpk2773avrg63927PLt3s/tmd709v3e2P9zc21/cu3n55BrCdyvwZL67IvEzZ/ynm70fPd2dw/gz0/0Wq3UxWf8DG/+Vm0srS1540N+1mYsF8WSPhZ64A+Lx7dH6xyfn5x/3KhzxPnu4v75ZLJYPtusTl4boKq1e7N7d31ydb0/Oo8M+jwyYCXTT/vZ8+IeXzz7YPlngHzC9f9wyffc3vUxHkX5UDxYnxg94I26XyWm4X3rg6+WXHPD0V1+it/gOftwS/f5veomY0K6vue8Mn9JRKTuRHV7lRPFbdtLfzhK95u448YWIay/vXRyu9yiOvvk/vXn54ttlE/O4bmie1/sv/Hl58eH2FRbdxu/VWMw7bnzpXo+PLq9frLYDpLq4KPpFEWef9Pado9y5WAhrn3uG3vQs/tDdE6/5OR6Gt0kzj6Yc76N47vVWifZ8+/J1icXGe7u84u2JhvC6n5ZR/3Rp9QYLdnvgoytuTEBR/ejZ/kV3JXx0xQUnvpCffne1OHU/ukIKL9bqR+8+7Kvz4Kv9zb+/3F77jRDSUhNR05KSERbbvLl6x7Vdl/PR9tqdNlgZOraoKk2CSCbkoGM0SblK01xbGeqYLERpUlRyLobceLS9/pB7W770+wvWEtiXl2gZLYcSmtQSqwUqfLLlWoi7kcxCiNMjmQRjPJIZxhyj5ZqatJaJgo25tBJqysmvDhrCWCTWUHIMjbjQWKJpshi0aCRIMpaaVbgOKbRCNGSsBDU01pBqh1OLzZrGIkS2w1hbzpajSuXqmSGMTaNKSKFE5QqaMLYSVERTLq0RfhqbWW4xa5KsBLRGiy0mzSlZ8ajNaLWm0kK1EP0qhFFCsSytxprIeh0ll9qkZA2JidZRjBuACPC0xOHTMcSaoqVY1DIJwGOoTVOJNZacyRYelXuYihapRS1SUKTmKiFJaOTwjlE0xRJiUiFipWNMKTapYjmS56BjbCXXmoJpJjt1TGoSNMdaJXJefkw1lBosNlOifzrmkHIillYjIWUdc65FMkElzT62IlBmaZlF00HHkpLmElIwK5ygHkurjbhV5M6jQccaQ7IWGhYZiZVjrVnNCmhVEnDHFixYK7EUI7ioYysx1NoklgpT19Gk+oTUJo2kaMgwSlYl/sbJdu6gSqEKscZqaYijRMuxSJYc1TIFLZYUY1NCVUMcg9YarNbacmp8wQ1KKXNTVUlxiKMG5jtXq6U2pSBbLlXNopBqF8coMapWlo1oehxj4iGHxiIFCkySEg8tMUc+STFW1dKyRhIe4phqaTGk1IISRY9jDla1xFjNmnpB0diKimrRykiKpNqsttBir1FSsZCrMss++NJaTVpUoFPGUqGWZi1prIluaw0pxChVNTWmowXVkCRbCWpGQdaiVolctszgmmmJrULphSYsqgarElk7xmY1RGhFixLG9Lu7hDhhtBJLpCC10IR5rKkmChqbLcdUIrlCaQyqliFb1ZILBdmydqIlVJrGYDkHZUML+zGNGsWkFCN2z+GFkS9bCsFaJRqfRmUqk7HeRC3TGGOMYqkkTYSu0xhLCaFJMiu5gEeEjYTUzCQrqCe1EqWlZJGAeBpTETWuZ6vWSKsZk8F1Qgup1UyjmT1sVpl8ZbQ51aoSs7IB6SXXDBWWxn4C0yLKqMxaqT5jReEzJVkpmSZKViuEZ3NqJPeMhZWUVixyI9yQxiqSYm0pSqqVGlW1xVoqxO+f1JSilpSLaiIBaKwlNZMmCc4BXrWlHGqroVYyINLYJIZlP5F4l8amYgScS2iEytNIFzkrXK0WMG05ppCb1dpaALFWoNgcSy2BTloLllo1SSEXb9NSFTWy00hiSaNJSbHWah1NCzWkbITHU6QL09JUUmmihLgTbsaUs5aQjHEZk5BqYecoU2HJckukTjRrXpCzFBNjFdQPalj2q7HQjLbn+9+fenwo+/XFYff8mqun0Bwwl16eb990cCzlb4lI/xgF5wHeqkdbLgw87Fyip9FKKZERtEg+huUxlRLVatZWtA5lbFqrxSASc1YX6ffPz29b+YRru26v2fIXKAt/2PzDld98tuE///NTvLEXz4f9xWH4dBmL1/2Ymwv/F8PqBubmweXFzWF7QbTqZvf0kvsmhRyAP7G5D3dbovNrU+FnNPXJ/vmL7Ulb+kc3jNcF23y5ov67f9y9uvrFr2TM//iLX7z61Zfg8Ltf/p//8aW38Ltffn/pL5mEEzVvuYqua3poo9vrR9trlDxXn3+7w0F/Qk7Yz+v0L/evPbh8cXV5s+8XeFL2RjLFj1aZcSr/TarMHhVYvBaCp+/PqjPXSg5WZrvH5rdFvmHx/KDOHMaQYk6GGEiRdDC01hoaXCWaF4TSrLWYW60prhryD+nMcayWVVooaGs/QWVOY2mqSZBDKH/wZURJLK0g+DtPDciShJpROodrYpZFTLUMeZSoEgMCuUZO4Y3Ic00iUlqLiQILTWLLJIXFIY9aXFkI2UrlfYxWcs41ZCW9MI9JWpWU0HDIEMpjqgk9tqQQNVMjx6IthkY/gU6zpYxUTtJIY85jSa1ZLNESqgoFZoqstFyLI15TzmK5CmbDkMfaao5Zcy2VvL08NkXahCLWHM+WW44ohWaov3lsrcbYqiCjOOk9IrCQRUXRiCmIARGI8oUqkUdLRcXvQG1JQMtySU00Sm5kXOXRikpiUpk0b7RoKxXFRsl94xP03JSyRNI1KZAIPYpARxTEYi1WjVabOR6hWUHk50DWeh6bYaBwr2oivTOPDfVGK6ZKrF4DpbMi3nLzRWgSgiFxvS0mrEgpuWkpNXq3NWSkraFCe0EpxbC4rHFdLIsQUJolxkIye0b7l9QUfQ3jKI/JiptGIpWExTymCEpNSy2ZY8FjLCkk7qWtmSTSPKq5EdtKv6k2jxpalprRakhtzWNAIFvOTZT81DxKkRoYu+RmaB9WawHpViBGVAfTgApXIlNGgSSNEjCWugKDjSfoTtqKK1aqTuC5sZyoXtowaYpkS42CFBVbIGNgVgoimzzWpLWSgIii6RZDIiOyq6LsvGaplKr+icQWsaVSIKsVFbnmyvW1oE5Bi0Uj2qumxLVSI2ujuaAlVUEPL9HV7txyaxG1O0fBZmcLZ7cglNHnkLOQ4BjHqDVpsiy5qGEwqJZYWq0t1EqCNZcCN6wYrRoj5oBoCVlSq00yp9tG06pFRGLC8sIEU3M9NTUW1804zalli2qNo3tjicUKRNli5HjzmGFBJaMxk2usY0rWcpYYQyYFVMeYW1bDoJfmxqOW1iRVqYb+jR2Ma0FY5kqCJ6ZzgkVI1sBhAx0FD0DCVIGbYJ9rjZYl5Bi60d9SwXpQMW0kjo7VFcQWqvbk2rEYpo+JsmlwLISCbpqKGIZxGHNKpA0XSSl5jVRjkYL91/N4R9wDIcUUG3Q8hDHGhkEEfXky5qg11RJzK9AlokSDVikZTs3KhjFkrTnFGmBxuCvEUrNsDeInYRcjFp6lpeQSBxmtZRTsXJgFrmkbLUaoIxUTJSotY2tJIvZQFVGuWYPoWIOcpUKag4zVJEPspjGRqyxjTQVuqSXhEaHlKglFPVWLUZsOMpaSSrAackoZ40zGokVLacUstEJeODm57a+uQH+vuvOmAg33jiFU7LTcU4xlTAUXlGHosR8xq4qwqlUrwvq/9eeu1v+X0p/j36r+TOY6EUCPKClBjJ+uQPtXOJ3Rv+86nQWt1kIpOLLIQv+pXmcZU8xRzP1wGS9QGyO/IZCiFRyrnEkbc4gWNcWYcvFN9h/4nSP8yaLAXXEL/gQlWvBHRq2tJWSg88NiqLABX2zjCAe8TqqpKiIlcvHPWKUioUvKIiXA63CzaAu5aCrNq+A3j/BWCSguMlatEtRyrGocXZKxxpgE/UaCkbJJiVXXi0PA/0FJymh/CUGHFS9jzRpiIwKQtXDd6lhzU22aNKWYuVMTn40ioXHFRvES91aqBvUMf8HRaqhZLZTml3qOqAQxS8zBQuYyT9Tq2HJGU8SnRImpJYl45QonLShpEWOgBI2cExDcQTWr4tMsOVOniYkkRR2reL9kbCGFLOhPSTnPQ4mF1gRro3AyS8amiZMMDRd2dFGFctMqJ2aEQ0IIr1i0WOKAS+4CDl26YpXE1LiEdGwJLbdh4qSa/KtUWTfc0jEWRzAHVZFipdUuOltGtSpJpJRYex1DnElpVooYLZcolgvO/Yx7lhI02ZC0JryYlOCkTrkUbTlm/6rG2jihkiR1Amu1CvI4xMIZOqS0aLDFfWU+COI1DX2/CAdoqFMEAwDJ3bzdxt35OddkeLCoYiEank0LeN0oiI1tgsKQek9oaUR+oqRl9fCGRm0WTFrXK1AV8EWyTysrY6Ko4YQUOOhIAUspxKREnSbxjFZJ7gUWHzduWNTZVNELwM9wfXsMIVUOglHSLIaoQSWHXmAtBRzLWJHUQFHPeFMFc5ICTBs2jabk2wrLMzKE2KI4/RnqYWME9O8IB+V3MIwl8B1DQcg1B6/BIlgITRvIWsjcd0xJgi+EiEt/wcYamxGHvPY6EhU9Tko18XbEclATqzk1XwbD1ouZUVps7E6TmkKDdnPs6AjmVmKrYWJ4lVQlWyuiibM8PukhmAZNuXFYipIAIhH6rt5OM7RGSVkS1o7TQE1EIVSLRd/BzTAyiA5VzFuvoxZKLTj4OcME5QgbFa6mRIehv1rh2pGQkS9ea9l/ICQUxaLxOsr61RSTceKQ7WASc0U/XfvCAmJXEz3xyWhVPfiYBMe2U38hNlUrNoN2Si45molXXHovoeTQshnqrPeVQSykHLUlX9KWcT83KSljUoJPwkRHh24cgaYAH7SmXEO0vkWSWCPcww+ZVG84llRz5RBbLuLTE7FbOdSnEccIzIv3ueHbyU6oTVUzNoP/IIfXCYwYfbwZ7haYoCoGAVykRB+6ENZD70+GpQAzDVUjcaCY1YdVcQg0Q4626jurmrTIvAtxNTCsLWuK+KcwT72kWs6FH08hQOBn44ag4a+un3+vOvWmfs72U+J1HgRKftqujPwqDVG7ymnGNsY+nYZ3QDwC/d8a+n89Df1wediev66k/0J++U/kh/wiLv8t7qU/SW78W0gVgduuajvpfCdBhB/KFUko/q62RxT4O3p7DNpUInppLJ6K8pP83m9L3rA1dQPHLmHrUPCa4bUK+K+SVngSbpoUCp5oZClnvSNu5uaabq7cWUfyBkofSneNnHFF1BLH5YeqmnMrfFaxhNQQ9QnPqRTc2MQEc3ZTosIHiFEWuJ/zsdV0eOsglkiLZ7H8YKUfnclCvkSyqCUUPzuL5wvPbWqJc7x4TiyTqJEQk+SUuCIuhUwBT9WQnFrVjCqJTqujJvJFyJNI2EiKNxL1x1rBl4e7TdHhU8DZyK1tiB+VVmrIlUwMnHo1J0NgNH4XBD9gTQj9WMkkIdcAj2cOVVWKpycUUggikePMqeI4NsnRGsfepWcjoMNiPDWTHoAOGX0e0ykUbp8YI6qOsCJ+kDqNKUtJLWS0fPfWEs0XvPfJ086Ia0AlyfAam4c1jPPmOZgQpPDYsSLyYkO1wgMsRmJHxa/ucXBL1iJZAA2lERdxjqhkRkaC+4wbVp4HiM3j4maF5BCMNm4EyCO/k4Y7HSMhOLlJ4+fTMtpCIxwgmiNHpnMStEx+I61WxbHsqHk0xUOtSWtxVzb+YGHVUQFwdjdNFbxYAvePg5eGHGuOBSyw8HLTSMg7UFAyKozmwE96uSMfnY34f7XItYojsf9ipVQ2D5+g9CdPEmmRG0fQr7KpqxAlgFhtkqQ1TJ/ChYsjOQTWai2lqNBoxI+H5ysv4QJhoWtIUFhkiksjMmBoerGnVZRSoPzS/Fh9xGmMxmn43g2SZHlYlCzWGj5ivPoZ91r1VKBMJpoIxi0ai47BUtCk5KgpFwvwA3gN2xK7lXsDyJAghiGxlOo+UnRTyUJODhGZMAb0PVymOGzZfqxfJPyVWnGV0CCi2mqKzUrX90rW6GQfF21YsOBabMU956hXpRWJmZ8v08RPWYw1xZJSDgxCE9pwZQIa8aaUFk+rJrIcStGIxeJ1clD8wpBZt7ijGNELxu36X3V80f1CJmTkHoEs7prBY4EeWUkICZBnbZyucuOeWA1js+4urhiGORC1oH3XPrHvoBuO/bsCTz4TZjfcKnodNxUrOSB4EPgqGz72lotTISVVLbn7WCuXorjb2VmVppK4OgODwnmVe4vVzRATrdhMrFzt1pXg2ycKAXUxLsIDmQAQng6fHhxRmlBiVSE86sREpK4WSV2rJjrWXGVPWGZuKCWCNx6HqtpNHrKhWsXNX7qXANuUKyvIJzGfDVL4Yi6VPKBeJ5MAKPiOMlkp3dtAZlolotLNK9JeWiUDi5Q+6kj0wFwkp8wNt9rIDYIQzDP6WNNa2KQELpu7VfovBAquHiv8JMtYhSBd8Z/3y8bYQY3wo8dMCnWwiQJJQfgyWIqikSBcMcVWGsTDr4hQwVfgBlh2LuNxZEncSUqC6F898+WnaHFvmAtvleVda0nuR1B8j2RZcu8NAalAiiHxIriipcKtKqRG9UCZaoIlNYsW+2Vjo+9cdZMcsYxDqMBFDX7uDoPieYxBI4ZZCHXM+H6yEjfUoY7ZyAPzMHKq5krLX9ZkueKU35pD83PSca5IvV4b+jnJOFdEaNaG9Gdk9VxhSq4NxZ/TECnKa0Pp5zSUTxrKP6chDn2uGJWf01A9aaj+nIY4lLpi1F7Pn7qS1zOnroLnTV3pmjR1UvDLf/ryKr7xQeof5OMHtwV8UN74oPYP2vGD24K3ZGRhNv2l8rEebb/ZcYzZTyR837E07LsT8xPDDfBnnmUMWsbCz/wGTYHo/MZPNsLPYkpkQJH1iJXZTza20YTaBNYTxf1XgcnZiaVpjRFVdHP2OPSfjEVbXD4wov7e/vqq0rGQJyypf6XLK1u+MqnkO3g/66vYv9LYEncericqb09Y/vqCKzaYoE++2j67/PbkBM+Dy2vOUm+f7fnJStbzeBjlyTUdYBj33y7lB519MThM/P236JAYx8vlCh08Aq+fHz9eoXN7wc1yjPCDS05Mgiixs9ufou6X3uwvfuKvsb52IuaD7WH3itF8Hym5+d9Heer36Mc93zraUbAe7wx5nYHTU/THEfNyuYlguYLgR5Nv96282B/6gaK+tByG+fh6/3x/cXtKC/vVzEIlxcvgT3d+/vxPmYXljNS9f76X7B7/8yl666z8rU1I3zgtFDbkn2VCSGK9d++f790LkkYEy4/eEevk/OeRx59/Nj45PLv37u6bez4j90Iau0P675Q+jr/h7HnDd/9x7dJvL1qO/358dfeeo++7fgERRmL0nYs8Hm2vd1wj9ofNP1zvzjZTyH/8wZ4dEe++/7Lzst/7VZf+S8+eCP3O/oLu3ru+/pjfoeY98Mff7K7Pzi+/5XIkkm6urzmZfFfFv70fk0Pi/Ijy6p+0QKgNnyKuE5xjHhnBodBaxG8VR+KohUMFkvAbkN6HWcjRD44qkRFY4uIYyomLO0mNTdwqnMfoCWg88CP2vSTlitMol8jNlXms6qmbpXGoA1+Rn98jbZDzQJyowYXKLfijSk1+BI6gqydkYWHX5tEpVbwZbkOmkj00XAsRqu5yKMbpG48RWe0hL3IKPWjbaiKBzM0STBrqkdJa3U7EHETQY5lnctt4SlwdTb3acmw9s4DbMRsYFKwgcTs04obE4KnCeSq+dU+MW/rko3GVIC6RWv2+yDGQqoE7Bu8C17hSv3J/7Jhwk+Dmw1Z3nxfnTGL+C9wf+OR8URD+dc/ZBNcv3nYrBdocPOFtEvjPqMxB9txKYWVUxQXl0VsENVdU1JE058LxQFwYR0VOOTtH7jEOHr/vqStyJM62jN/Os19PFDltkSwAPLS48k4VOYv4TnBaxMIFmCiLi7ZWgycOkyVccYSiMK6KHKnlOE45tVbZqq9ftPMnKnLrOi0qzsEvsPgJYkv/REUuMry/BUXux17V8fcitbJTF4dRo4ZU7soOdroncT293l9xSMZ1Wn7A/rz/gP3Fxe7p4faOrIf7V7tnXaL4T8j/8f8DTDP49A==').then(json => {\n",
       "   const obj = Core.parse(json);\n",
       "   Core.draw('root_plot_1779222316447', 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_1779222316447();\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
}
