{
 "cells": [
  {
   "cell_type": "markdown",
   "id": "625899d6",
   "metadata": {},
   "source": [
    "# rs_bernsteinCorrection\n",
    "Example of the BernsteinCorrection utility in RooStats.\n",
    "\n",
    "The idea is that one has a distribution coming either from data or Monte Carlo\n",
    "(called \"reality\" in the macro) and a nominal model that is not sufficiently\n",
    "flexible to take into account the real distribution.  One wants to take into\n",
    "account the systematic associated with this imperfect modeling by augmenting\n",
    "the nominal model with some correction term (in this case a polynomial).\n",
    "The BernsteinCorrection utility will import into your workspace a corrected model\n",
    "given by nominal(x) * poly_N(x), where poly_N is an n-th order polynomial in\n",
    "the Bernstein basis.  The degree N of the polynomial is chosen by specifying the tolerance\n",
    "one has in adding an extra term to the polynomial.\n",
    "The Bernstein basis is nice because it only has positive-definite terms\n",
    "and works well with PDFs.\n",
    "Finally, the macro makes a plot of:\n",
    " - the data (drawn from 'reality'),\n",
    " - the best fit of the nominal model (blue)\n",
    " - and the best fit corrected model.\n",
    "\n",
    "\n",
    "\n",
    "\n",
    "**Author:**  Artem Busorgin, Kyle Cranmer (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:36 PM.</small></i>"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 1,
   "id": "3ef95cb2",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:36:09.064003Z",
     "iopub.status.busy": "2026-05-19T20:36:09.063888Z",
     "iopub.status.idle": "2026-05-19T20:36:10.046105Z",
     "shell.execute_reply": "2026-05-19T20:36:10.045543Z"
    }
   },
   "outputs": [],
   "source": [
    "import sys\n",
    "import ROOT"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "05bbd7ff",
   "metadata": {},
   "source": [
    "set range of observable"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 2,
   "id": "db44d86c",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:36:10.048583Z",
     "iopub.status.busy": "2026-05-19T20:36:10.048439Z",
     "iopub.status.idle": "2026-05-19T20:36:10.158898Z",
     "shell.execute_reply": "2026-05-19T20:36:10.158457Z"
    }
   },
   "outputs": [],
   "source": [
    "lowRange = -1\n",
    "highRange = 5"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "4245fff6",
   "metadata": {},
   "source": [
    "make a RooRealVar for the observable"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 3,
   "id": "9a762a91",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:36:10.167481Z",
     "iopub.status.busy": "2026-05-19T20:36:10.167348Z",
     "iopub.status.idle": "2026-05-19T20:36:10.319686Z",
     "shell.execute_reply": "2026-05-19T20:36:10.318583Z"
    }
   },
   "outputs": [],
   "source": [
    "x = ROOT.RooRealVar(\"x\", \"x\", lowRange, highRange)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "32ad5d88",
   "metadata": {},
   "source": [
    "true model"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 4,
   "id": "109890ad",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:36:10.321178Z",
     "iopub.status.busy": "2026-05-19T20:36:10.321052Z",
     "iopub.status.idle": "2026-05-19T20:36:10.552023Z",
     "shell.execute_reply": "2026-05-19T20:36:10.551368Z"
    }
   },
   "outputs": [],
   "source": [
    "narrow = ROOT.RooGaussian(\"narrow\", \"\", x, ROOT.RooFit.RooConst(0.0), ROOT.RooFit.RooConst(0.8))\n",
    "wide = ROOT.RooGaussian(\"wide\", \"\", x, ROOT.RooFit.RooConst(0.0), ROOT.RooFit.RooConst(2.0))\n",
    "reality = ROOT.RooAddPdf(\"reality\", \"\", [narrow, wide], ROOT.RooFit.RooConst(0.8))\n",
    "\n",
    "data = reality.generate(x, 1000)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "a704b06f",
   "metadata": {},
   "source": [
    "nominal model"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 5,
   "id": "6a6d1d75",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:36:10.553764Z",
     "iopub.status.busy": "2026-05-19T20:36:10.553630Z",
     "iopub.status.idle": "2026-05-19T20:36:10.760000Z",
     "shell.execute_reply": "2026-05-19T20:36:10.754845Z"
    }
   },
   "outputs": [
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "[#0] WARNING:InputArguments -- The parameter 'sigma' with range [0, 10] of the RooGaussian 'nominal' exceeds the safe range of (0, inf). Advise to limit its range.\n"
     ]
    },
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "[#1] INFO:ObjectHandling -- RooWorkspace::import(myWorksspace) importing dataset realityData\n",
      "[#1] INFO:ObjectHandling -- RooWorkSpace::import(myWorksspace) changing name of dataset from  realityData to data\n",
      "[#1] INFO:ObjectHandling -- RooWorkspace::import(myWorksspace) importing RooRealVar::x\n",
      "[#1] INFO:ObjectHandling -- RooWorkspace::import(myWorksspace) importing RooGaussian::nominal\n",
      "[#1] INFO:ObjectHandling -- RooWorkspace::import(myWorksspace) importing RooConstVar::0\n",
      "[#1] INFO:ObjectHandling -- RooWorkspace::import(myWorksspace) importing RooRealVar::sigma\n"
     ]
    }
   ],
   "source": [
    "sigma = ROOT.RooRealVar(\"sigma\", \"\", 1.0, 0, 10)\n",
    "nominal = ROOT.RooGaussian(\"nominal\", \"\", x, ROOT.RooFit.RooConst(0.0), sigma)\n",
    "\n",
    "wks = ROOT.RooWorkspace(\"myWorksspace\")\n",
    "\n",
    "wks.Import(data, Rename=\"data\")\n",
    "wks.Import(nominal)\n",
    "\n",
    "if ROOT.TClass.GetClass(\"ROOT::Minuit2::Minuit2Minimizer\"):\n",
    "    # use Minuit2 if ROOT was built with support for it:\n",
    "    ROOT.Math.MinimizerOptions.SetDefaultMinimizer(\"Minuit2\")"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "e8b48c87",
   "metadata": {},
   "source": [
    "The tolerance sets the probability to add an unnecessary term.\n",
    "lower tolerance will add fewer terms, while higher tolerance\n",
    "will add more terms and provide a more flexible function."
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 6,
   "id": "ce3b96d8",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:36:10.766412Z",
     "iopub.status.busy": "2026-05-19T20:36:10.766282Z",
     "iopub.status.idle": "2026-05-19T20:36:10.994614Z",
     "shell.execute_reply": "2026-05-19T20:36:10.993905Z"
    }
   },
   "outputs": [
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "[#0] PROGRESS:Eval -- BernsteinCorrection::ImportCorrectedPdf -  Doing initial Fit with nominal model \n",
      "[#1] INFO:Fitting -- RooAbsPdf::fitTo(nominal_over_nominal_Int[x]) 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 6.28407 ms\n",
      "[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_nominal_over_nominal_Int[x]_data) Summation contains a RooNLLVar, using its error level\n",
      "[#1] INFO:Minimization -- [fitFCN] No discrete parameters, performing continuous minimization only\n",
      "[#1] INFO:Fitting -- RooAbsPdf::fitTo(corrected_over_corrected_Int[x]) fixing normalization set for coefficient determination to observables in data\n",
      "[#1] INFO:Fitting -- Creation of NLL object took 159.78 μs\n",
      "[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_corrected_over_corrected_Int[x]_data) Summation contains a RooNLLVar, using its error level\n",
      "[#1] INFO:Minimization -- [fitFCN] No discrete parameters, performing continuous minimization only\n",
      "[#1] INFO:NumericIntegration -- RooRealIntegral::init(corrected_Int[x]) using numeric integrator RooIntegrator1D to calculate Int(x)\n",
      "[#1] INFO:Fitting -- RooAbsPdf::fitTo(corrected_over_corrected_Int[x]) fixing normalization set for coefficient determination to observables in data\n",
      "[#1] INFO:Fitting -- Creation of NLL object took 120.911 μs\n",
      "[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_corrected_over_corrected_Int[x]_data) Summation contains a RooNLLVar, using its error level\n",
      "[#1] INFO:Minimization -- [fitFCN] No discrete parameters, performing continuous minimization only\n",
      "[#1] INFO:NumericIntegration -- RooRealIntegral::init(corrected_Int[x]) using numeric integrator RooIntegrator1D to calculate Int(x)\n",
      "[#1] INFO:Fitting -- RooAbsPdf::fitTo(corrected_over_corrected_Int[x]) fixing normalization set for coefficient determination to observables in data\n",
      "[#1] INFO:Fitting -- Creation of NLL object took 114.48 μs\n",
      "[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_corrected_over_corrected_Int[x]_data) Summation contains a RooNLLVar, using its error level\n",
      "[#1] INFO:Minimization -- [fitFCN] No discrete parameters, performing continuous minimization only\n",
      "[#1] INFO:NumericIntegration -- RooRealIntegral::init(corrected_Int[x]) using numeric integrator RooIntegrator1D to calculate Int(x)\n",
      "[#1] INFO:Fitting -- RooAbsPdf::fitTo(corrected_over_corrected_Int[x]) fixing normalization set for coefficient determination to observables in data\n",
      "[#1] INFO:Fitting -- Creation of NLL object took 122.171 μs\n",
      "[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_corrected_over_corrected_Int[x]_data) Summation contains a RooNLLVar, using its error level\n",
      "[#1] INFO:Minimization -- [fitFCN] No discrete parameters, performing continuous minimization only\n",
      "[#1] INFO:NumericIntegration -- RooRealIntegral::init(corrected_Int[x]) using numeric integrator RooIntegrator1D to calculate Int(x)\n",
      "[#1] INFO:Fitting -- RooAbsPdf::fitTo(corrected_over_corrected_Int[x]) fixing normalization set for coefficient determination to observables in data\n",
      "[#1] INFO:Fitting -- Creation of NLL object took 124.45 μs\n",
      "[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_corrected_over_corrected_Int[x]_data) Summation contains a RooNLLVar, using its error level\n",
      "[#1] INFO:Minimization -- [fitFCN] No discrete parameters, performing continuous minimization only\n",
      "[#1] INFO:NumericIntegration -- RooRealIntegral::init(corrected_Int[x]) using numeric integrator RooIntegrator1D to calculate Int(x)\n",
      "[#1] INFO:Fitting -- RooAbsPdf::fitTo(corrected_over_corrected_Int[x]) fixing normalization set for coefficient determination to observables in data\n",
      "[#1] INFO:Fitting -- Creation of NLL object took 124.37 μs\n",
      "[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_corrected_over_corrected_Int[x]_data) Summation contains a RooNLLVar, using its error level\n",
      "[#1] INFO:Minimization -- [fitFCN] No discrete parameters, performing continuous minimization only\n",
      "[#1] INFO:NumericIntegration -- RooRealIntegral::init(corrected_Int[x]) using numeric integrator RooIntegrator1D to calculate Int(x)\n",
      "[#1] INFO:ObjectHandling -- RooWorkspace::import(myWorksspace) importing RooEffProd::corrected\n",
      "[#1] INFO:ObjectHandling -- RooWorkspace::import(myWorksspace) importing RooBernstein::poly\n",
      "[#1] INFO:ObjectHandling -- RooWorkspace::import(myWorksspace) importing RooRealVar::c_0\n",
      "[#1] INFO:ObjectHandling -- RooWorkspace::import(myWorksspace) importing RooRealVar::c_1\n",
      "[#1] INFO:ObjectHandling -- RooWorkspace::import(myWorksspace) importing RooRealVar::c_2\n",
      "[#1] INFO:ObjectHandling -- RooWorkspace::import(myWorksspace) importing RooRealVar::c_3\n",
      "[#1] INFO:ObjectHandling -- RooWorkspace::import(myWorksspace) importing RooRealVar::c_4\n",
      "[#1] INFO:ObjectHandling -- RooWorkspace::import(myWorksspace) importing RooRealVar::c_5\n",
      "[#1] INFO:ObjectHandling -- RooWorkspace::import(myWorksspace) importing RooRealVar::c_6\n",
      "[#1] INFO:Eval -- ------ Begin Bernstein Correction Log --------\n",
      "degree = 1 -log L(0) = 1216.78 -log L(1) = 1208.89 q = 15.7692 P(chi^2_1 > q) = 7.1557e-05\n",
      "degree = 2 -log L(1) = 1208.89 -log L(2) = 1203.21 q = 11.3692 P(chi^2_1 > q) = 0.000746732\n",
      "degree = 3 -log L(2) = 1203.21 -log L(3) = 1198.85 q = 8.72171 P(chi^2_1 > q) = 0.00314444\n",
      "degree = 4 -log L(3) = 1198.85 -log L(4) = 1190.16 q = 17.3777 P(chi^2_1 > q) = 3.06393e-05\n",
      "degree = 5 -log L(4) = 1190.16 -log L(5) = 1183.56 q = 13.1965 P(chi^2_1 > q) = 0.00028048\n",
      "degree = 6 -log L(5) = 1183.56 -log L(6) = 1182.57 q = 1.98352 P(chi^2_1 > q) = 0.15902\n",
      "------ End Bernstein Correction Log --------\n",
      "Correction based on Bernstein Poly of degree  6\n"
     ]
    },
    {
     "data": {
      "text/plain": [
       "<cppyy.gbl.RooPlot object at 0x55f953065190>"
      ]
     },
     "execution_count": 6,
     "metadata": {},
     "output_type": "execute_result"
    }
   ],
   "source": [
    "tolerance = 0.05\n",
    "bernsteinCorrection = ROOT.RooStats.BernsteinCorrection(tolerance)\n",
    "degree = bernsteinCorrection.ImportCorrectedPdf(wks, \"nominal\", \"x\", \"data\")\n",
    "\n",
    "if degree < 0:\n",
    "    ROOT.Error(\"rs_bernsteinCorrection\", \"Bernstein correction failed !\")\n",
    "    sys.exit()\n",
    "\n",
    "print(\"Correction based on Bernstein Poly of degree \", degree)\n",
    "\n",
    "frame = x.frame()\n",
    "data.plotOn(frame)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "0cb4cb5e",
   "metadata": {},
   "source": [
    "plot the best fit nominal model in blue"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 7,
   "id": "4b9fc419",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:36:10.998291Z",
     "iopub.status.busy": "2026-05-19T20:36:10.998159Z",
     "iopub.status.idle": "2026-05-19T20:36:11.173650Z",
     "shell.execute_reply": "2026-05-19T20:36:11.172979Z"
    }
   },
   "outputs": [
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "[#1] INFO:Fitting -- RooAbsPdf::fitTo(nominal_over_nominal_Int[x]) fixing normalization set for coefficient determination to observables in data\n",
      "[#1] INFO:Fitting -- Creation of NLL object took 161.541 μs\n",
      "[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_nominal_over_nominal_Int[x]_realityData) Summation contains a RooNLLVar, using its error level\n",
      "[#1] INFO:Minimization -- [fitFCN] No discrete parameters, performing continuous minimization only\n",
      "Minuit2Minimizer: Minimize with max-calls 500 convergence for edm < 1 strategy 1\n",
      "Minuit2Minimizer : Valid minimum - status = 0\n",
      "FVAL  = 1216.7779341626624\n",
      "Edm   = 4.16187387343553302e-07\n",
      "Nfcn  = 19\n",
      "sigma\t  = 1.18138\t +/-  0.0315451\t(limited)\n"
     ]
    },
    {
     "data": {
      "text/plain": [
       "<cppyy.gbl.RooPlot object at 0x55f953065190>"
      ]
     },
     "execution_count": 7,
     "metadata": {},
     "output_type": "execute_result"
    }
   ],
   "source": [
    "nominal.fitTo(data, PrintLevel=0)\n",
    "nominal.plotOn(frame)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "8b00d744",
   "metadata": {},
   "source": [
    "plot the best fit corrected model in red"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 8,
   "id": "384118fa",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:36:11.175115Z",
     "iopub.status.busy": "2026-05-19T20:36:11.174991Z",
     "iopub.status.idle": "2026-05-19T20:36:11.322792Z",
     "shell.execute_reply": "2026-05-19T20:36:11.322396Z"
    }
   },
   "outputs": [],
   "source": [
    "corrected = wks[\"corrected\"]\n",
    "if not corrected:\n",
    "    sys.exit()"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "727cc7f8",
   "metadata": {},
   "source": [
    "fit corrected model"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 9,
   "id": "c2e8ea79",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:36:11.336409Z",
     "iopub.status.busy": "2026-05-19T20:36:11.336266Z",
     "iopub.status.idle": "2026-05-19T20:36:11.487956Z",
     "shell.execute_reply": "2026-05-19T20:36:11.468725Z"
    }
   },
   "outputs": [
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "[#1] INFO:Fitting -- RooAbsPdf::fitTo(corrected_over_corrected_Int[x]) fixing normalization set for coefficient determination to observables in data\n",
      "[#1] INFO:Fitting -- Creation of NLL object took 272.222 μs\n",
      "[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_corrected_over_corrected_Int[x]_realityData) Summation contains a RooNLLVar, using its error level\n",
      "[#1] INFO:Minimization -- [fitFCN] No discrete parameters, performing continuous minimization only\n",
      "Minuit2Minimizer: Minimize with max-calls 3500 convergence for edm < 1 strategy 1\n",
      "[#1] INFO:NumericIntegration -- RooRealIntegral::init(corrected_Int[x]) using numeric integrator RooIntegrator1D to calculate Int(x)\n",
      "Minuit2Minimizer : Valid minimum - status = 0\n",
      "FVAL  = 1182.56856160934331\n",
      "Edm   = 0.000128899179831313442\n",
      "Nfcn  = 185\n",
      "c_1\t  = 3.18336\t +/-  0.783721\t(limited)\n",
      "c_2\t  = 4.8037e-06\t +/-  3.11006\t(limited)\n",
      "c_3\t  = 1.18953e-06\t +/-  1.51773\t(limited)\n",
      "c_4\t  = 0.966776\t +/-  2.14197\t(limited)\n",
      "c_5\t  = 0.215629\t +/-  67.5047\t(limited)\n",
      "c_6\t  = 10.4446\t +/-  20.5657\t(limited)\n",
      "sigma\t  = 1.26666\t +/-  0.209614\t(limited)\n"
     ]
    },
    {
     "data": {
      "text/plain": [
       "<cppyy.gbl.RooPlot object at 0x55f953065190>"
      ]
     },
     "execution_count": 9,
     "metadata": {},
     "output_type": "execute_result"
    },
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "[#1] INFO:NumericIntegration -- RooRealIntegral::init(corrected_Int[x]) using numeric integrator RooIntegrator1D to calculate Int(x)\n"
     ]
    }
   ],
   "source": [
    "corrected.fitTo(data, PrintLevel=0)\n",
    "corrected.plotOn(frame, LineColor=\"r\")"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "a82681c0",
   "metadata": {},
   "source": [
    "plot the correction term (* norm constant) in dashed green\n",
    "should make norm constant just be 1, not depend on binning of data"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 10,
   "id": "bccae9bc",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:36:11.498832Z",
     "iopub.status.busy": "2026-05-19T20:36:11.498675Z",
     "iopub.status.idle": "2026-05-19T20:36:11.609068Z",
     "shell.execute_reply": "2026-05-19T20:36:11.608622Z"
    }
   },
   "outputs": [],
   "source": [
    "poly = wks[\"poly\"]\n",
    "if poly:\n",
    "    poly.plotOn(frame, LineColor=\"g\", LineStyle=\"--\")"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "0d0d8a05",
   "metadata": {},
   "source": [
    "this is a switch to check the sampling distribution\n",
    "of -2 log LR for two comparisons:\n",
    "the first is for n-1 vs. n degree polynomial corrections\n",
    "the second is for n vs. n+1 degree polynomial corrections\n",
    "Here we choose n to be the one chosen by the tolerance\n",
    "criterion above, eg. n = \"degree\" in the code.\n",
    "Setting this to true is takes about 10 min."
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 11,
   "id": "bb823f44",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:36:11.611286Z",
     "iopub.status.busy": "2026-05-19T20:36:11.611144Z",
     "iopub.status.idle": "2026-05-19T20:36:12.317493Z",
     "shell.execute_reply": "2026-05-19T20:36:12.317140Z"
    }
   },
   "outputs": [
    {
     "name": "stderr",
     "output_type": "stream",
     "text": [
      "Info in <TCanvas::Print>: png file rs_bernsteinCorrection.png has been created\n"
     ]
    },
    {
     "data": {
      "text/html": [
       "\n",
       "\n",
       "<div id=\"root_plot_1779222972308\" style=\"width: 700px; height: 500px; position: relative\">\n",
       "</div>\n",
       "\n",
       "</div>\n",
       "<script>\n",
       "   function process_root_plot_1779222972308() {\n",
       "      function execCode(Core) {\n",
       "         Core.settings.HandleKeys = false;\n",
       "         \n",
       "Core.unzipJSON(44028,'WkwI7jUA/KsAeAHtneuTG8eR4P8VBmI/3EWU+iqz3t2xHyhKOvlOD4Ykr8jVOhQQiSFxGg64M6BEecP/+8UvqxuDGQxF2pa91q5IAuysR1fWKytflfiP1bf7n15uLtYvNqtx9dWD9cUP66uvN999ebF+efV8t1+51dnvL7b//mrzuw9Wo3ers/e3+6v+9Pl3/2/zZE/6imKfv9xvdxcz8H+3F09XY3Crs8Obxv+4q603NRBC1BCjW519sr3YPNid7y5Xo8zgl/ufzjfX4Nfbp/vnHfxoe34+FwZZwKWw91Z9c7b/dH35bHuxGv1AyhfbZ89vJb2/2+93L24W+2r38mbCo7MtSKhbnT2+frzfH3nxo6v9ek8rrVHmBnS/Q1T+6HL9YnMbb9JudfxQ7maHDkWX5OWVVD8eGV75/u7y6ebyy+0f59E7Svx093TT5/WR9P8fz/8/0j62j+f/H+1397+7erh9vTn/fq6x390FrsbcMl2/USEWG7LrKjcTVuN7PeH3N2oxm7+/rnMMLs38y40K/SX/cl3lZsKhmaUr+92juTfWlVOQmfReYsilSvZNWmJ5Xtd//H0fp7k+IGgu4Gp8jxeoSM0x+lxKy9nWyfnux88+eDAP+jHw6PcvLYNSj4+evz6kfnx4uv/d1Y0X3f/u6sa77n93dV3t/ndX1zV///qF7QZG+Kfrx9cv1q97h37/0+Hxq+eb/Xo1BuvY8+38dP/q5ebJ/ov1frvrvfjs1YvvNpf9+avtk+9fXz/+1B8/2T2bEz/ZPbtO+2PPfbh++nC9vWD3uNXZg8vd1dXz9XZ+4QF8uJsp0fGyZnF1+HpFf7p7uj3bbp6uxrP1+dXGrc7+9+X26eub4E/X4P3vrh7sdpdH5T98ut2vv2O/7y9f8YKPtq83T2/0e3n1w8vti+1++8Pm6oTkfbK9gqIu1HYG15eXq/GbP7jV7uWehz+51dmHrzdPrlbjxavzc7c6+6yT5ycCtfxquwePGfrs1YuH6/PNfr/QSAbss83r/WnqB7/78uEn9x+vxtU/LY9udfbB7tV355v3X52dLRP2xWa/3l4wWnNfH11t/7j5/dWS//gmaLlfbNbnq1Fp3LI7LHHQmgrL9+vtxdPdj1/tXrKPjuHHx/BMr64LfLyBOs/L4Mdlpz94vpq384P1fn8yzPf3+36Q0b9H72/2P242FzOlvgHZUH50uXvx1e7lapSBtfPo6XoPITTg8QJwkt3vgPzJrb7/dPfD5vOX639/dVgR33+xYVxuJp59vH32/BO6MB9KtjrX+yfPl8H9/svnux8//GFzsf9yv96/ujoswu/vv9rvWAaHkp9uLl69v77sMIvk/hOW2aHG2Reb9dPPL85/Wmqcfb3dP9+92h+vyGWVfry+mtfYknJc6ptbZ/UvxhFAMd/IEXy9+c52/fbi2ZvYAlbGg/P11dW8JSjX+ZDjhJcs1ZUfNSU3fyYZvfPOT2qpPIXRW66fouVpSlM61PFTnkuSXuayPNexNqcirsap8VyjU8mT0GB082cSGaWpk6pOik6iY8nO/k0SRvHq5s8kcRQtbv5MkkZJwc2fSfIopbr5M0kZ1Uc3fyapo2pz82eSNmpUN38m9VZYGvX9pHIT1FFacNLESa6ThlFqdlIBw6Sxt5uTk+QnTaPUQI6TlCYFqxnJ6Ccto6TkJKqT4Ceto4RkHax50jZKSU59BpMp+FGCOms6hynIKKGAgZPsp6CjBO8kVceLQ+jt1NrrMlTa20lhCmkUKc7qxzSFPIqPYGAohzKKb84qBJ1CHUWjjav6NoV2PeopTdGPUryz7jaZooxMseEc0hR1VBVn4wPIUNWOZPBTjKMNMbOmcYpptDEFDWlTzCOvadG1PMXSm2GgYpliHa1FMTSn2EaJxVknRKfkDQlbYlOSUWJ2tmanpDYOMxDGVpfnOJa4PKcxyfKcRy3LcxnD8lhZ2PNr2uhdZAGnKfvRu1z7M1unSX/WkU7ZtpoyO0hC6zlsIlaHVU8AtXQg04b0RrLtpDC3Yq3bJpUp07wBWqZC+wb4OBUQAJDqpwIGBqQ8lWUPM78FDCzH16mAAUBNUzEEUnJZprLs5FCmsjQvYSpt7Bg6P1U/dgR5ljFZN3jWsdow8MwW7tTFT5VVacMA0Lcv7fmpsk3yArB1O1Z+qmxcGwYAtq2NtvNTmwlXqADSKZJEgE68pEG+WrAcyRkgdiAqQKdiLCo/tU7Gmr25WKkipFd7joXnZs8anJ/E+1FEXExGFMSLrdCoLuokXjtRaEyvTOLDWIsTWzVkR+uHpOwCYBpVqwvFpTiJz6Mk7wobtE3iIRnZykrmVWxPdSE4VUCGxGVGVSYx2prYQzqJQMbEBXEKpDaqEKPcJoG28sxH0yQQV2aTgee9kjrGfT0L1NWLa82J95NIGSFw9EfqJFBX1l1uTqDobZTaaFgmUT9KhGZFJxn6Dh3r41YBdbRqvrjaJtEwarB2/CSQ1po4FhzbXKCtEpyRMqVwtq2uUjvOUNejI0C9jE3Ykk0mZUogif0zqYc2BTd/JmVOhP7YZ1KfOr2FrE/KpNRgJMBP6svI3iziCm+uo4bW0ZIyqW8jtSDWaVKmhFVA5UlFxiKuV55U1OjyXHtSJoX2qDypRGr2ypNKomavPCnzAUWXQu1JpXRsl0broXPWajN8l2Y57iBHHelJbUrAtrfLlPS+Sg2T2pRYZzXESTUeemstMylzb0FaoeLXo6pl6W+vXA8dpnJY9jC0zccpLNuYxcS5xU6eT/EO981s+V6nwH5e8g0+4gkM7rua8kmmwMaeixtYDxyEgdC2zm4kmaK/Zoz8FH3nCagOpIeaQNd8ClDkOOFEmaJnNDr1on/RNrPRfutPZDMDcxqnMEWWjtdD/6L4DpPvdYowS6DICWqw0tbcvSjBIB8NfzE8KEpvxA45Khu0kHo/RTFS3/slkHqqkM5BQwU/ReWc6V3SeUo6zzhF2KG5SxxpkU07d8m6zPqYu9TzoSRqXe75cB9zlxgCrUuXevG29MlKBw9IA1Y4iIF91qYYlkPPqgY79fp8TzHYwXuoaCfvoR7DYZM2JZbkvEZ4S2pymOoOX089GKV2Pfk9P/b+9D5OiSXZ2zGUU8ug3GdtSq0skL3L+t47OKVmXQfyOmXPFNipOGXfuXTqZE+fbdKm7OmxTdqUfefW4RhYhX2HGYoZ8tV32AxX6/I8bVPmTCG/T9uUjYTZCuxNsgr7CrQuZY4W22H9dRLsdXQRBDlZbBHOhft4LN2SY/Ehz1w7VI+afXfMneNUYTX24VBbCMzalFmSS/+oB8Wa+8cqyZCsvsNmeD745+nPRrN6/2yojUu3HTbD/XjrCyRrW7pH23DpvXcG9f05lwzXQlQOfVDgSDLcOXvMnvtw9GdbGjyWhRB2tKcyE8J52KfCqptRTDIVFt0MUjkd+muZ19NP98rM3hgKU5kZHFqywr13PbMuNNCoQ4UGzuOaZKq+UwKblKnORNCIRPWHOSdnRqfTnjqvxhnfqS6r0ShNZS3Oa8/aOAw2PGanyL2JeRFSmIIcpTNZoaRRwvmV14SQnJkc9e5WztGFIHudKudoX/ZWuB4ga+NAi2BybQHSPMB8OvTRmSpH6ELXeS1EktXWiXGFRs4glWdOeEZpXo3zKEx1Xo1LS32vAoESMuPcEMz2TLs6UizGuW9si4rIOA8R26KGPkh93Uz1mkCCcbDzApx62Zlqzwd1DYt80Atf02xryESUpWo8llFqNIrNe+y90Sh233xTjSYpzeskLgSbfsZrdMAuXqPTX3NAh47FAzq98IIOL0oLNgbMI9SP6KkiLs7zYh1Jh2Vk701gNFPn2qXGeUCnmjpG84BPNYGR5aIKSYbRNQhGM32sJkPO50HN/UDr0zbVLkf2SZuqCZIz4axdkOzYTjX3470vhanmo0OW9k2cvH5rH5+l1yZRHqp2kXLJK7MqZEGozKffglLpC3sZlNIXdqfYtfSFvXS0dGboMAoFyWaRsGvpC3sZCATNnsmKqn1h9+Ov1j5tS9HaZ60fGrXanPVTolbjgWZ0kDhtbfXJReacaaaNUe1bf2mjb/1DG7b15ybafPZYtZlSLyMwE+q5/Zk5sGOqzkR6BjpRXGZzptGH+Zxp9PKawwFkC3qm0faiNpPoGbim0IxbOyLRrOAGkZ7PYq9Tm8l0b6XNdHruSVvo9DzVbabTfTDbzDMsUOcY5pluM6GeM2ftmq0pkELcYOD7Z2oQS9Zx/0xNbQb7Gdo0jikgVkugqi3ziGhcpqZ5TN4k8IjQX8bYENQl8tI6xuoQPGOemrYxZpe8k1inFvwYIZ9OkkwtyBijM2Ti1IKOMTgjT2VqIYxRHeJ3RpMQx+gdGresUwtpDM2lhi5saiGPobrcVYQtlDEUh4hcaKGOITkE4EILDbVJzk5Knlr0YwguFyelTS3KGBRpHtm4RR2DiXxS09RiGLU5ZOFaphbjqNWUBA21RkJtVRLawqmhScuuFCeNFgoagtKcNFpAWHWVlSxTi21UdbXzdc30aM5EvzI1iGNzNaGTmhqksbqanUqYWkJ74Gp1KmlqKY6SXUN7VaeWTCBEEaYyNehick1niC1v6gnlnfB3qFw6ZGJzyx3K8FWulRliybiGnI+qBp0BugiVNrUcgATFhYGGjK1aA1PPDXNu7mBXXLRcOkgn69RypVvi6SUgClmHasLAgurToSpVKVMrYiDqMANR+ThUQh20ERLRGURL7NBeoAJvUMbqxHQNgGgXur7QwNJB2iUXYd4JWg5GG8oICBppatUbiKbAQOkgaJCL7tqUHh0MHQSNODUUcuTSfUAUsk6U7gPmDoIGYOkgaISpQSApDBqAzcAAGmFqrWOF7kl0aghRFe3jDPaxQt8kMjWoJNpp0ADsMxhAA7DPYAANdHC2oCSABqCtKNOuGmhLSiLru00NkT46QYvl0aWhtANmwVdglFEQjAVGq+4koo6w/DAKGujs1BfKo8FEI73AqcN1gVHLOIkNjTjli8EJdAzmVIHsOPWJfFOnSgIfYPFjm49iA8VAM2OQqyPK9K4rES/BQHCx3GggqARy04hKB0wMzIAZRNA4osBDn7iA1UDQsNw2os1jVFAxqre6ph0EFBqCnDXGVBUkMX80hsxYXVO3SmPI4HV1hukT3C46SKgTaJsSb4ZB1NjdbuZo4ILKCMyhX4ZMX/5GI2keAaxhtnFSaT/Ascww7SOEYVoyqwQwGsauU6zgE9DFzjD4BITDrhdFP+mDWQRMT9o1uQXCK6hnK6MRKmQZww0UXDzknfaBwQcCT/ud3ouHxNM+MO1D5GmfE4D2Y+Cw6HpMxiNGjpIOg09MHDQCPS+0HzPHkFTGnPZj4ZDqsLVfOcK6tYrxiI0DrsOMR/Icf2Z+y+CThMPR1KKZ9pNykHaY9mee2GxczEeKY6I+BzrtwxcXJxXBgfZTHjkjDbb2ypjBh3xrr46sHYPpf2ojJi70vqYqzx4NYodpL8uIzpr8SHvZDIYdpr0cRsYGC12kvRxHTIsG015OI20ZTHs5m6LbYPqXy2hae/gIa6+OGP+ATR2f29iw/gHTXvFYlTrM+HIseBCG3tEgBwMCjiXQovHMDJFRRBLiyDHSE2iT44FNU5uZCcVzQKBlJUGtVQgLzZJgzWIWoFkSrFnsduxRb8Y28RwUimmRBGaWo4Klawk0y2FhFkdKsLY4LliMWD7R83sODFYjJYRmOTJYftcJ6OSxWlKC8efYYIFZAohxcMQZD6wFnqODJWYlQIzDA+7MEhigRQdnCWBqSrh5PARMTR8yj4dRejtEGFMsFqBux8g8pkbLOUhs5VECTDlKIEuWAKYcJp3OYGokAYZgnlsSzAhklKA6gcyaGYitbQlYXRZLkCWYhQa7w/IOTCmcKWaIqA56ijXIVjMwthUzPdC16iCoZhBiIGrBUCpmEDIDenEQWCxC0lhvxUFgsQlJo1vFQWAxCwl2HGCzEqHGAaEF7poqg0HYzBB0CRYPGMkFfLJZbDANYbWkPAQY2xAmYYPBD+OQLeaKkVGwDmFkoT8QXAxEamu7ugp+6FKg+JhzaR/t3rzUIbhmIoLCN+8guFiJzOYMTPto96D4TVxhvNCn2EYQB8HFSqRQeGBrHwsja1odBFfQ70HhgRkflCq2TYKD4AoHDKatFlw3pOmotmuig+CKHTAs+OgguMIBY5sIph84jRpY/4gEwHnEhmkw7XPA2JbKDoIrHDCR9jKGWZHQfReYUwiuRD+qbbDiILiCSR4KDkz72ORtv1UHARbkCNtu1UGABUkimnTgEu0jS9juaw4CLEgTxqCYiCmCPGF6Bu8S44FEYbA6CLIgUxgcXKR9pAqDk4MgS0JGhbXIDoIsCS2pMWAOgiwJ7S2cSHMRfBI2WxhHMdOqJPADDg4CLZjsDY4OAi0J/MjPDgItCfyAC8ZskQR+wM0F8Mngh7TiHfRaMvgBq4NcSwY/4OgC+GTwA04OYi0Z/ICLC4xH7uPHGsJIKRn8vFNYY/DJ4AcsDsotGfyAg1PDB/yAMZtjagU/4Owg24JKxmCkIWDwIx9OHJj5bY41AxGX0ucXXxcBn5LmfPYgcJ5hOBJg5EXqw3UBs/6AMe7iSPT+qzN8yf7JHPJWH53v1vugK7c6N++xlNzqh9X4TdPomiL6ZdcU2lJd0+Za8K5h4g3qWgiuhehaSK6F7FqARlXXQnMteteiuBbVtRhci9G1iPCIjwi0rboWm2vJu5bEtaSupeBwN23YwFN2LUETq2upuZa9a1lcy+paDq5lHE2Sazm7lqGl1bXcXCvetSKuFXWtBNdKdK0gimbXSnGtQIeba9W7VsW1qq7V4FqNrlWE1OwatNTod3OtedewQTd1rQXXcHBpybWGCFtca9B9s6ujh4ESe9OEeCwQphn12GIgwR62wUN8PQ47HrLrobUeAuvhqz2k1UNPvUnFUFIP+fRwyR7C6eETPCTSwwl7iKOHInrYAw8t9BBAD9XzMLMeeuchch7e1UPevDkuQcg8HIGHhHnolof99FAsz5LzcJuexeahSp5l5rvzDjVYVR7FjIfm+EQN08rBVHrOeI/6xcNTek54DzfpOdo92hbPme7RtHgkEs+J7lGueI5yDyvpOcM9AoLn8PaoUDw8pOfo9nDqnjPbw4N7HAg8HKSHpfYc1R4G0nNGe1hHD+/pYRg9XKC3o92YBo5lD7PkOZA9fKPnKPaNGrCNnpPYwzB65hytBUwDXzg5MOfdq4k5RxHBec8XEj9zbowiWgYxDrH7QzDnxh8Kc94dA5hz85dCNcCpypc5UtAGc47cz/HIF20w58YQdo8s5tzYQWHOEdY5u/iiBnOOEM55hCKDNphzY/4QqDlY+KIGc24OW8bnIR5zLPBFDebclIHCnJt/nHm8CHNuDm6IsRBqvqjBnJuHnUmh3TOMOTf2zTzkjG8zodAYtu5Cwpwbu2bSmfFpJndJqX/405/+5P5W3pzp57w5+32En7k/Mnts4qZ9eUHBxbO5p9w7vODWRQjAo8sOt+85XN8AEe5dnF4B+XR9+f3m8uhKSU84euWccLgm8dXm9f7+xTPcrnFABeyZfvCMgeWfb59drNCKdPjo/WR/tMObPZuH8fr19tQ3/P5+f5903Lefbn/YXm13F1erMQktknP0wk/W322Way60Z3BvIdKCwZ+fnV1t7P4JdHZOPKAdDO/tk+8/2Vw84/aMHzxOzjYHS1XrC+7mt6vtzxd380ORpXkckx//anpo0/kX9PBffzU9PEzQnzmH768vj24Pvb++XBaFuWBzH4ttev7wy74nPrhc/9ivXXT485f76yseHZhveXRgvujx+cv9B93fvt8gwwmeTWTb6POX+5kk0InPX+4/shtXc9GPtvO9gBMPegpY4tPtnmtjC/zVbnduDvQk9CsqD3YX+92ry6v57sL9/YzOLYp5f79nExuR+hlaoH8mMWCv0P1+t4hOAnEjYjDow4unH15e7uYLW+xsA604TX306uLJTBbIBDyiYoDzFJLLlZe5MP0HnAuz7wGP5vuTzbPNxdPjezVg11OPKCwvuk5c2l5u2vGKhSwcCjJ+81J0q7OPuQuxubpFw+fUL1+un3AXwNo+XH876sPh7tucBo6HcjexORRdkpeit5q2crf7fZR4fano4+0VC/IYH5J434xO9nR7Kbc03EdnKbqkzgVvYUOpT7cX2xevXvzr5nJ3fdWDjBs3Eo3E91svDy83Z5vL//3JdemefjRwPeG4m2B6nHrdz576webs49WYPJN2SPl6NZabKY9W2AGOijyeEx6ujxffw/WNtUXjh6Trli3p9Cbmw/XTG11n7B6un57e7Hy4fnrH5c6H66cs9kfXwzOnPL6Rwtk431Wiwe2T7+ebSg/XL/sVykcz0TgkPF5hnlydffnkcrO5+Gj9xKgP6EHWjoYfkH1wtGxJOp6PpdbR/qEI4PXuIWVZQH1Z9TKXL6CiqzxEu0hDotGMVmdcjMQgnHbUvl7uO1HyYwCxfuwvty8/2DzZvlifXx0uFRlJnrkdPbAMR72zAre6Z2nH/WORWOJRBw/w4ajqTAoMSL/xe6h1IIFWp3enHUrTn+X9dKcTn93vLi42l1/QPUqy1ey1V6vxGy4S3bvHl+g9Qei/F+/VBc734j3JQOopEOci15AlHaerv1n8kFf9nFf9vejvxXvRSlrbf6OvP0BcNuunm0tOa7szZcN2gD7a7j9aFk2aF43dxGIWDzm2lmzUnqzPrTKz/3922wsSF0bgwfrlMfjV9sWBnSy1NqnRKMbvXqyfbXjRgcA/WF88Pd98/Xx79f3m8ov1xbP5jnJPf3/3ek7rs9dTDZOja5r/st2dby+W1Pn6Yi/6YHv55Pw2tZ+zuF8K0kcH4CNY7g9fv3x0zPYsiY+PEx/fVXJJvFGSgp+uX3+wfWZ36FmEn1/un+8erF9sLtcz9TkV1x6un/5i9++YsTfev3u4fvomeQ1JihE6oliAC/GZc4/G79bJeCqFnZLrO2j1bxfx35OhJM+flGpBl7iyi/nv5cHnFkVL0lpyZIK4qJ/uKE5yHLLGgjKlasgsvZvX+CXnAY2K16wlFNbJ0a1+SW2QRG6pivLyOnc1Nj+0GFqrueaAEuX2jX/1dchKiRICqkErsEQMUB+G5smtoUU7opbgAKvxvTCUoNy3TylJsdybgQHy0FJrrdUYY8D28nNxAqoOvjKYAQUTmN4MGqBBBrUCpqKz7AOemoeUqYzG/jprNb53O+s6JgAxBd6ToWpOwedckjRh8OegAHP2yQwv2V0S8K2lUkqpXrWcxhxIaYCwFvUxZBRqh9cTguCuhTK/vock0Fxqi5ol5FDVpvYQg2DwnrFtrZWSY9XEy6/DCgzpRrYtmjlMAVEvbNl6rxqS2KT/TF6PTeCHZqPvvY+xdCVFD1Tgh3iScyPkwZ2o3giCcAe2h6AIdzV8CJJwV9tz0IT3GK/TqAlM8Rw1IdYhpaStRp9sxf2ZIRRo4LcQCocINX92CIVbkvy71z8EXfj2VtiFDr9r4IVf7qL+zM78Xe7ps+reyCeYOPwmTgEK8ctyCsSjsX12HJCGHcbBdrK9joUMetHhWag8Za7+Ss7qo99/9oCR6PGOfnbUPpYP3jRmyQJnzAvujEA9377+NqWzxqGRpCEhdJZ9Na7u3/tit3t4vtvf253d+7fV638zhv5oyLva4a9kz/58JflnTzbniIkYfFZnb1ByzxruN4WUmofgtWlXj/r8GuBXrRYH/++2aPQxja3OHlmgHVvWjyy6jrFkvQSxYM4+2l5ezQrOT9bLE8G3dBFuX2w+2F69PF8fBa5B2jrIa4yZWQKuA8l8unv6yfq7Gf4ZTf27TdNPt6fJ4phc3fu3/3Xvf9zzg8/3/uc/1LwddAIH/cJsBfkZc8bRvB1mrc/fEh3p7z9pbzA+vNuk/fH2pP1DTRGDbLqJX/kU3bKXHBlTvG3/Dy/2l1s0Tyywr65evfhx3tU8Ljuc5yVylj3PGZ+uX6MRvuPwm3XFnJfdhPLZ7vLFoohk2c72jh5v6uzL3pjRmxvh+zAdEM3v1Ez57vxTDyl1raxcQkxZa+9vL9A1fnh5+TlBr0AN+PMfNpdn57sf54G5f3kJdbvpkpPx1TGXHDtpCNxkjjk3OYFf3uz9ko3y9nP+i90O5TyF7zrlzNvkENHr+beXm/X5dv/TB+v9miq29lfjilfsnl2uX3DIH5X59uX5bv/tt3YeHh35bJu/4sBn+H/+wIeru20V/+zlbntBMMj5RPu5ifLdd+o9P7Ti+Ba+a+K7NPsOfGfLzZabLDdabrTcYLnBctVyxXLFcr3lenFmlfeD8CTF+YHAJoM254eQnB8iGQQ+GbiINyQyMhkERhkKGZWMSgYeLrxUeL/QlNCqgICAi4CWgKGArIC30AWhNzLgnEEfhe4KPRcGQbG368Clm4F4OgM3WwbupAw4zA14vQ3c1xjwwRu4DTFwcWEgJMpAOKWBYDQDQXaGRkZrLmAHCvQ80PNAzwM9D/Q80PNAzwM9D/Q80PNAzwM9D/Q80PNAzyM9j/Q80vNIzyM9j/Q80vNIzyM9j/Q80vNIzyM9j/Q80vNIz+PQisW3e/wOi4ULSLgzOsLk4ENoDwS3ssTqcDwsjss2eDwSlgXXweZCjyJTrZZSkftIXMLpcU4czkS4uds/nG9ccTgBOlF7qk5cdMUFl5266KpTF1xwaumkAol9xP56y+GbgEs9rT9bbANL4+n673Wud2Jj8ksQ3QPdOMQOPDkSDqfHckbcoPwfovx5++R0bxH8XbxJuB8+er599vwvqIcy6e3VvgkDjmktxtBCxjkqDjnF5DW0FGNsmdUXU/W5tRY0BErUHLOo+lBjUluQ2RfJWXGVCy4NPmcvsamElGzZntQ4aeSkxK1W06CpFXz/SwsxJ5eGGHNSjSnFkH1zaZAcShENNaVUs0tD9Spao28hxaouDSGFFiTlWJNdax1aEZSA6EmDUiWFknIKmqpPpbk8SCiSAy54VbPc9Y6sRUIKUqRmT/eLpOAlVyk1WJVbqMfh5oCdjPnJCIYht5xQr9aqUZQNn3MIMeei5uV5xzTFEGMsrapKr3IrAUpUqsdtM/pQjSBVlRpqayEamWMcs4bgYyiQs9vwrcUTBs0plRa9b8nu8g2pSiw+puxzNiLcYogZF1BNOFGe9AximjSn1mLw+KW+/R06lFxDzT7lmLUoR4wWjVnVt5Lt3GgSQ2hSs6pyDfU2YjLkoMXHXDPel3ZclNaChNIUNTEnU5MqdEWqOV6evvQWHndUudXKSbMnVd7el5N36HCC+tuavaOVW1VOSvwlCaeVTgbxjjKnLd2RdHvK/5IEO0ZOZ+COdx0dOHdlnyCzlL+V0WPSPn43Gv9NHDTnFlvJuHNzl3womBdy9BF/WyhcTkmqtCJS8axP0IkW1OeCYYYqNScNMbeUK+7reZAcg0af2UBa76py0srJS281m4fA8ZE9lq/so8tDCtXHHIrmXLiZMmj2GkONOSWsL3motYr3KfhQk4jLA6YKaTlinanJlQFjVUxQdAh+hvLSRKilxtxcGVRqrCFVHwtu+yevyNDmkkL03hfoUeaMSDmXHJNww+QE85MBOxn0kxGMQ6UTvqSUaxVe6nNO2JKCKB7xd8xTimgkVVsyz+d0muArWsuUo5dsx6z6UKSFioWvwhPmVlIJjLFvxlueJNxaQHGIkpP6pEEyhxfMqtfqzegVuVEz4Gnuo4TCKMFv3uoc7Gmr0nwqoeEi/vZ3hKHFJFGStJrBVAdtDaqfGpZU3iFZC4eq0pl4+lIdoMY11aC1BgqkGnKUFFOrkVvNQw6lppY1xJQMr5N3nqBxq8rtRk7hkwpv68jtV9yB982O3a5wMlbvUOAWUne84lYJGSSWWmLMcIFcw79d4mSG3l7l5B1/QZWTd7w94aSVtyecrJ2TKreb/YNb/SDEur41bqsLWY3FrX7Q1fjN7Uonbzlt5+SFb084aeV2FU6db39Cr2lq7NV4t4L520PQ+m97oPqYbhJYt/r2Ynf54l/W54Tr9ga+2F6sz9/fXswuK2irSf9y++zF2jRK324WFV1X2X17uf7xoLVDb//t7nL7jLd8bTHKTVV2fHXsWE81qz/+W0m0v7it7fydFXAPXl3+8EYDZdIUju1tF30tfIuq9JvXfzjWwj283PETJ9vdBWo4co50btnUZm/WuuEN8bNeUG/XuqHW+xmtG/4rP6t1y7PWTQafuWdnD+4905+1aBq4at9qGrhs33xly0yWmSwzWqYlB/tWy1TLFMu01xJQbyDi3eAr6jdAIM1o33girPEQSUfllkjPpGfSiV098FBJb6Rzp4zXyoDmDNUNreKng+qNdPuQzoUrsBY6IINp30jnmtvAbTZ6zRlBIJyBEDYDwWcGNEqDKd1IJ4jKQOSTgQAlA0mZr0J6IZ2A7wNquIZg6JWvyknD8Q6kpKOGoreB3obB5FLSM+mZ9EI6D5X0Rnqr6NvQXiDnCnoNRYxW0k2rQV40fQfpiXTTh3QFB3wP6ZX0Fl3iL/3ly/v8du1bdtc69W+kDHDkLYbafLFIT7cT6hA15BJ89iVyXxFXsQpjX1vOqEXbUIuUpLG2zAU49QjoGqLPUjPDKoMKTJGEpna/V4YaQkoaY8mxECxpiMgFXhK/KaM6tNhariV7mEqnYYihpqiVy54WPHioJdWknTFDiTgoWoLSkPDRJQ6pRtEcizR6pnGosZXia1TBhq5p8CmkUnPydU6QFuBiS03cHE6D5takwfdzCZSE6ktNyiVL1AhpUM3JK91Ti6E+cDuzthQRdwytppFxkorHE2jk0qIiGqXMaKFmDTkVn4QPCd63nNVrrZXYUWhyvZQUW4UFQgFcg2+hxJC4ba1DaiVJVQlomEjwQVIqEmsomQBEQ9QYQ65wxsRi9kOpPqCUiCVYWKABZUhKPldfuOzZhtCazwg+0e7g1iEzeDP376QMjYVRTcFmoWlgbH0NORR67yQPodErqSWY6iINWUKQEIKWij4ELj7mWEtMicvJEgefavLVpxq4oBkGrTWnEIO0xP1aHZJWWOucsoVXkKFUYeYl4CPoRAaffCoRHZXapdsh4PmGd2Up3MoeMvqVwu3PmnN1dWigY6q1HD0JIaB8rOrh2qMrQ5GSkdFE0em4MogWRJbC3lBkq5RTKbEIo1gR+nxAZqiSSo2Cyg6RMYumwKK3fZtCaAG1lFaixQ3ZF421JB9CSEYhassIVwXfP6hI9V60hBBb4CIt+6I0lEmhVK71okkFlJRitQBRQ2nK9WrlAisBxIbkc/LRN7ax0UbWhoZAO4HlMngPRahFcdk0WsuOSDHn1Ozu9JAqI5VriMkusw7Rp1pbLKGi1ISARy1eJUcuMvMO31JWlLwsB7O9ZBQANbZQY+j6sph8EU1eakqZMgVNqTbNjbngSDEROTDUISXhCElZYo0aSpZCSAkc+krUlrK0Krli8YmqPkkrIbeihfeEHGmn5lSkREsRdr8mjY3lYcdaK/wcVFBpXNHG1iTSUkXPIUyiHyAv6ptE9aVYW5JrrDm3WqRyo5pzMmgpKdScpCbwEUW9GLTmKIRAYP9p5SYwt9EzP/gw+JrZ9DEUVj8npi/BSwks29xSwIrlsyhrRFk1bDvMYxI0No0VnbANkY+0lGPShiqW9nxIBX1MbhHVgzWoLecSMxe60bxTSmOqWVJCe83VY5J8UM05JImN0AhcZ8slJQ0+xyrJWykJpTCviNhcYKeUBC+stBJ8EUPVN40pFx+j1ERcBTw9C1qbGGuLVcqclgXtMD7GKREDg3Kx1ViyRGVZMareRx+ytMaleJTwlhZQfKu0GnJfVP7tSf8wJpubEtIJy3oiHZ0cUjelo18L2/5kd3m5ebLfPP3zGfeA08VvjPt7ds2Lb75+Y9w5LP/rMO6mr8OyVjFuEmYWxeNxQhqqNGxF2Tg8+LUIYa+BWFAEehpaCL6gbse2CJdYKgEuQlPJ2N1xJuCwkVB8QEiKQ2s1tQYPY4x8HkJsMRsXJNlpGZK01GqDOSMuYx2STxhiMdsSArINwUdshV5a8hj0hwbjmyIqaYIb+SHEjLGv1ZBaIAH1aW4Cn4grA44kgZMJ9h35zMPsR48d0JdmJUJonJ5etBZcCJAoQmvZuLhGQlQ03TXUNGNaSkvKGQ3zAepeU2zYaUvBI6MMUjJXVWDPEELyoDGnosaLE2wzDRobxuEgEV244lJRcspRq8eMDfPuU2ja0JqbaNRqROTQBKugfqixorWHobVYlUPxHF/YPqPFMxswmQi/vNcQIKRggoatTE0L8VyYlFSzz7UiBDhJA0pszYEbPkQniQOWlpIxW1uMmzBoLB4Wi2mCzeZKQ2y5cS+D2CVDaMg+yFA9xs2QYvWtZmEIxbWhcBwTRyZKIygRx2xNGugqkWTqEHxrCCVodGtxZYDbEG3GLIWAaQUTi5QYpRBsDXuNTyRqyphDjatWLcW3kAmrlIZUSsoltgTjOdvfa461SmipIU6XWrGcw2AR2SkOMVYvIWHRIW5OHCR6n2KDtbb4WAMsD0p8zAAZ55nsW46YeWCkzOGmSk6+eZ9yQs6XEln0mJdZxTrA/sVStGmKBCwdqqRmvJVnX6J2YMKa97kW2zxDkqTJJAA4WZQWGMI1EB6IOF2KjJti9B6pw3QVEgM+DrHl4FnWg4cGaAyqzSILDS1gS4fRhRtFUxJiyDmbSp3FMJQoISDDRkYS9UpSyfCEvpbu85SNjU8F6wvMPYIll4PgIImrM0SPHBJFMadQImjGBKKMibH/mszmlZtoYFEOUprHc6LEYuodIhLFkFKIqRD/aEBGgUL5HPGbwLksYR1vJmQkONYm6mPRmjzmIxRLiSBBMed+YQ0Jomlqgd0jhblEpmB5ZMUZINnBmGu35qWIJIqQESSXENhDkAiEjMru0Viqz8Tmwbss16App5Y84Yj8gFBmSogcW+X32IaYpOWA/Brg/kmRoBAyZsDUZ6HkWmOpQWFlKRKimv0p+yzsXz8g64WGTIj8TMe1Cn4BVbDDmS5OU2T9NE1VCRLoB9UWECZSM5aYFLQ4/MAty8zUeKglIgw3DiY2FoIBEMsqy5ffrMCWgKhVtDWpilQggfhOIFBKMtmJuEapYOLVnDMDRggu7sRlTRIa/fINEZB1J1qIJ4aI4wnlFEts0A1j/9+e9Ktlzk9Orl8nc/5yd/7T2/ny9zeXF1f7zfbiHuVvqdejhSA44tIXpn222PymXreNx57lY2SBp9/U678i9XocoumONWtCuwbrcTuhFU0tB9TjaJDTEJOvRbPUEKo511S0bhUerxGLLw+qKEIkBVReASZJ4fFD43RMJJRSG+q+xP+4qHA/tdRqauuID4vEVPCsSFWgxIV7wqhRggQOKFcGU/codB7+wZUheo0++Fo4kGHW0KmVFDUmj0BRONtS8ZyrNUQcYUJpKXJO4HdClaCaoq+izaeEGlUjnJtkKZXbySDW0PmHAptHAt4lHkdGzj861wgw6PE2SloZD5gsrLhSheikOOTgLYCgEqJU/G9MPYWDjRLzNg8B7hsFUU3cesa/MqLAbaJZYSIxZXBOVw5jvJHMhwPNdVXQSHCEUkIS7tsLCehQU0oeL6DCTNJgSLWlwNXg1D1nYFVQA+LTQmuweJojoYTjUGLmMAxctA9mBEIayDEX/B9t+cDVBY8CksCygxbTXWtpeJXAynYPGZOmYG29DyGr9bc7U1ZUp02JFJAzKmUkQPjYZrwPHuNBg8aIOSGy2xL+MAynwmajc87wJ5k2CZcaBtxW4RckCPFWMY5V3DfhuYgIEAYtGUkx+Ei4SjzXmUIT0oJtafU1WgwALdVc2LHFBGImRrUCMM0+lEjPzBdecY31JSY2Cx3BQ1Q9TFhjs9Fo8gGvIPTqJMDmtuYj2mQPGtEnbaZ/V5w8wxBbTDj2akS5b4y/T9h2KkpoSqAUTzVX9MjZLHzYJwT3H61mzPOIxzEQfxNrVhzgUTOqUpyGmbXoo2IfUhzomCWiN+QcU8SL2LybcvEFnjFp5idgBnYlM+9xJIosNyWwKLp/T6zaNJQQEbdC1kgg0IwwhTMwliGiuOKGRnyC1HKSZJSkoqOtKLBDSea5hu+yOQ5j8CoDUhimBK8ZvXUdVDRE4kw0ZBlXhxKqNpTRrEqHdI+HX0Faq4SoHZpXgko0TGQW6hSf6pghbqiPES01iCaYVZ8tYubQIJi1tlDMkKSIKAUblzaL6hmGVAShusDDcwFgiCnWoD4xsvwWxRCIetBKkWbxYTFGoX0WrDBYpwpVmEKPqzVB2c0ZN6SEz7lFuoaOIutXwY6GgsTHUpCmctbKb8PQSipFgkeg1TBU1VjBhYAAKARCSCUHgWhgfs5DCxUXR5h0pMA65MBulZa0dm1HhMokdCbE0Q3cGIkJ827BvuEClmKcB5MoclHIQyxQEhwo8aUMFXNVK94jKhVxEeMLLsQh1GqGJUyKJUeJGbtOdhFVRmADY4vBCfttCb8KHvyEBT1Rl5/082/Mkd9/9Lsvv7z/6Ydvdnb5p8vN2WqU9otr6r87/+yDd73N/nD9g8X4As277rn9ggGD8JMhBgA6hCowOkEt8G2PhJNwEQk+h5YDuiRK62rE5i3VS0oZzsECl/RAOJDV6KtWT7hu3k23ifKhhDYWmBwzc/L+OYurAvwJLcIkWRtzFiyUaWh81hwQSx4vWa310D0Rf/sCZsfRCejX7y6IYgkWXz5fP939eBRf6cHuknBl66fbV1f9huThos5houb78ntuMDMTxOu6O1Dtktmj1Cq43A7RdohSex1Ddo4W9snu4tmGO+gmk326vny2vbDQNFxaR0T78+6J3g7Jsd5v7ELjXesIxOdeWgetp28KgvCu3V9G4zho3aH3ZM6B/+aIf/TwaGKOZNLl3uUsk3YR9cV236/c9ruvrIjPzVXwcKcZbhN21JTeKf7pD2614xrsNxevzs+J2PyHk5jNvwUBO6w6Bv0kAOSjM4sNyww8vn683x9tl19ZaFniElHmBnS/Q1S22Ca33faIyHFrERzKLSHQevzFQ9EleXkl1Y/Do/HKY2JwKHgjXolRJ4sgpZ2SIH6pBZkkGAqXZ3EG8UgMJeBm0okfLNNRhU6UVmMZ/HFxkm9GAStcTqtG5ohHPxdYImDlChOOCR3vI4tEtcQIswhjSDOtQZtrpT+P97slEtb3qzFGuNh3jAJ2HCOMAFt1gGHEZaxkftLJgnwdvfwXjQIWc/mbRAFDF5oIb4+szwAtYb56kLDjGbOZXLJ7FLAakDwx0QQsWCviUy2BvL5fjYJiIOIRV6Ri+Di83qKA3bVSltq2kDzKffjJooqFi7XxznHA7goT9m6RwOJJlLC/RySwu/D9q2OBsVtPQ4FZONw5Flge2FF/ZgAwFvtvAcD+UwOAEQ/2EJrhiXT4twBgd3P+UI5bh+WbOSY79m4G0T4NFQrzz+Y6jv9l2wqGvu+p20fpjVP0FxeUWA5vjwTy1cfyESXvYmsZpJmtvVq/eEkU3A+2V/sPX+8vb8QBof7RYL7NdemdhvPnI33wipM7B3NoLxMbftnIXnTwv0hgL7tU0uN6GePWI0TZ4688rtc/1CQxoP8VQkT9FsXLZvEo/vsSIA/W08LjHW0j0v7uu+jdo3jZouxhvK6f0UIdMl6vRhzGuE8guXlfC306xPXK3AYiAnKuyAicDyeBWk5Cufxdo3vd0tqgACRkl/1i0V1H3C+oAZyFa/gA3AX5owlfwGJhQIWw14RLIVRHU19NxEcJyF0P+8N9Ek0Wmtp4hlSsdG5mXGNxLZq+0sMML7oZYztmTV9Z4jS37GMm69FB0zfXah4jS29nUR2GGeNQI2aom0pAFsFfqAS8pIEjPmKZjHfVgolJFjd/huJaDfYzSsC+OP/zlYB3c05vVIEOPmRQtxJdD2q783oEFrDTJKC/sSIw5GLRtaXYBacTrd9fohqdLynf+2d+noI/rJE3j8qS+Y8yIH0j4TtZQ/5FBuTTzfri3r1/vncPA4L90Od/59H4cv/03gebH+7NIyK+C7dv3DX/2OvjhubcoT6//jIDpP1q2Bxh8vOXN39f7K6fPeFIO0TQXX5A5+H6csPP9/3HbHuL9U83W7am/lYBL49+g1aP7+bGOTpdjzenf4NferxavzCz0t9E1jWyNJtzeD4Wcy2E9F+pNfhNzP0lf9YRPqPz37+JuX+X3940wsv2OFgNfwtW3X+76r+IuuhXLubaj+ua+IqchbDXZdnVmMrgpbaSuO3dZbGDmKsaiB+BAazikPiPJ+YeOIr/lDO9R4olSAlxZNXizh5Hil0C9En3EuCHAr58crl9SfghO0U/3j57fr599nz/YHdxsXmyv/5Nw4+2rzdPOzE5W59fbf70/wFVdpxE').then(json => {\n",
       "   const obj = Core.parse(json);\n",
       "   Core.draw('root_plot_1779222972308', 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_1779222972308();\n",
       "</script>\n"
      ],
      "text/plain": [
       "<IPython.core.display.HTML object>"
      ]
     },
     "metadata": {},
     "output_type": "display_data"
    }
   ],
   "source": [
    "checkSamplingDist = True\n",
    "numToyMC = 20  # increase this value for sensible results\n",
    "\n",
    "c1 = ROOT.TCanvas()\n",
    "if checkSamplingDist:\n",
    "    c1.Divide(1, 2)\n",
    "    c1.cd(1)\n",
    "\n",
    "frame.Draw()\n",
    "ROOT.gPad.Update()\n",
    "\n",
    "if checkSamplingDist:\n",
    "    # check sampling dist\n",
    "    ROOT.Math.MinimizerOptions.SetDefaultPrintLevel(-1)\n",
    "    samplingDist = ROOT.TH1F(\"samplingDist\", \"\", 20, 0, 10)\n",
    "    samplingDistExtra = ROOT.TH1F(\"samplingDistExtra\", \"\", 20, 0, 10)\n",
    "    bernsteinCorrection.CreateQSamplingDist(\n",
    "        wks, \"nominal\", \"x\", \"data\", samplingDist, samplingDistExtra, degree, numToyMC\n",
    "    )\n",
    "\n",
    "    c1.cd(2)\n",
    "    samplingDistExtra.SetLineColor(\"kRed\")\n",
    "    samplingDistExtra.Draw()\n",
    "    samplingDist.Draw(\"same\")\n",
    "\n",
    "c1.SaveAs(\"rs_bernsteinCorrection.png\")"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "ac9a574b",
   "metadata": {},
   "source": [
    "Draw all canvases "
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 12,
   "id": "2965f732",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:36:12.326263Z",
     "iopub.status.busy": "2026-05-19T20:36:12.326124Z",
     "iopub.status.idle": "2026-05-19T20:36:12.442857Z",
     "shell.execute_reply": "2026-05-19T20:36:12.442529Z"
    }
   },
   "outputs": [
    {
     "data": {
      "text/html": [
       "\n",
       "\n",
       "<div id=\"root_plot_1779222972440\" style=\"width: 700px; height: 500px; position: relative\">\n",
       "</div>\n",
       "\n",
       "</div>\n",
       "<script>\n",
       "   function process_root_plot_1779222972440() {\n",
       "      function execCode(Core) {\n",
       "         Core.settings.HandleKeys = false;\n",
       "         \n",
       "Core.unzipJSON(44028,'WkwI7jUA/KsAeAHtneuTG8eR4P8VBmI/3EWU+iqz3t2xHyhKOvlOD4Ykr8jVOhQQiSFxGg64M6BEecP/+8UvqxuDGQxF2pa91q5IAuysR1fWKytflfiP1bf7n15uLtYvNqtx9dWD9cUP66uvN999ebF+efV8t1+51dnvL7b//mrzuw9Wo3ers/e3+6v+9Pl3/2/zZE/6imKfv9xvdxcz8H+3F09XY3Crs8Obxv+4q603NRBC1BCjW519sr3YPNid7y5Xo8zgl/ufzjfX4Nfbp/vnHfxoe34+FwZZwKWw91Z9c7b/dH35bHuxGv1AyhfbZ89vJb2/2+93L24W+2r38mbCo7MtSKhbnT2+frzfH3nxo6v9ek8rrVHmBnS/Q1T+6HL9YnMbb9JudfxQ7maHDkWX5OWVVD8eGV75/u7y6ebyy+0f59E7Svx093TT5/WR9P8fz/8/0j62j+f/H+1397+7erh9vTn/fq6x390FrsbcMl2/USEWG7LrKjcTVuN7PeH3N2oxm7+/rnMMLs38y40K/SX/cl3lZsKhmaUr+92juTfWlVOQmfReYsilSvZNWmJ5Xtd//H0fp7k+IGgu4Gp8jxeoSM0x+lxKy9nWyfnux88+eDAP+jHw6PcvLYNSj4+evz6kfnx4uv/d1Y0X3f/u6sa77n93dV3t/ndX1zV///qF7QZG+Kfrx9cv1q97h37/0+Hxq+eb/Xo1BuvY8+38dP/q5ebJ/ov1frvrvfjs1YvvNpf9+avtk+9fXz/+1B8/2T2bEz/ZPbtO+2PPfbh++nC9vWD3uNXZg8vd1dXz9XZ+4QF8uJsp0fGyZnF1+HpFf7p7uj3bbp6uxrP1+dXGrc7+9+X26eub4E/X4P3vrh7sdpdH5T98ut2vv2O/7y9f8YKPtq83T2/0e3n1w8vti+1++8Pm6oTkfbK9gqIu1HYG15eXq/GbP7jV7uWehz+51dmHrzdPrlbjxavzc7c6+6yT5ycCtfxquwePGfrs1YuH6/PNfr/QSAbss83r/WnqB7/78uEn9x+vxtU/LY9udfbB7tV355v3X52dLRP2xWa/3l4wWnNfH11t/7j5/dWS//gmaLlfbNbnq1Fp3LI7LHHQmgrL9+vtxdPdj1/tXrKPjuHHx/BMr64LfLyBOs/L4Mdlpz94vpq384P1fn8yzPf3+36Q0b9H72/2P242FzOlvgHZUH50uXvx1e7lapSBtfPo6XoPITTg8QJwkt3vgPzJrb7/dPfD5vOX639/dVgR33+xYVxuJp59vH32/BO6MB9KtjrX+yfPl8H9/svnux8//GFzsf9yv96/ujoswu/vv9rvWAaHkp9uLl69v77sMIvk/hOW2aHG2Reb9dPPL85/Wmqcfb3dP9+92h+vyGWVfry+mtfYknJc6ptbZ/UvxhFAMd/IEXy9+c52/fbi2ZvYAlbGg/P11dW8JSjX+ZDjhJcs1ZUfNSU3fyYZvfPOT2qpPIXRW66fouVpSlM61PFTnkuSXuayPNexNqcirsap8VyjU8mT0GB082cSGaWpk6pOik6iY8nO/k0SRvHq5s8kcRQtbv5MkkZJwc2fSfIopbr5M0kZ1Uc3fyapo2pz82eSNmpUN38m9VZYGvX9pHIT1FFacNLESa6ThlFqdlIBw6Sxt5uTk+QnTaPUQI6TlCYFqxnJ6Ccto6TkJKqT4Ceto4RkHax50jZKSU59BpMp+FGCOms6hynIKKGAgZPsp6CjBO8kVceLQ+jt1NrrMlTa20lhCmkUKc7qxzSFPIqPYGAohzKKb84qBJ1CHUWjjav6NoV2PeopTdGPUryz7jaZooxMseEc0hR1VBVn4wPIUNWOZPBTjKMNMbOmcYpptDEFDWlTzCOvadG1PMXSm2GgYpliHa1FMTSn2EaJxVknRKfkDQlbYlOSUWJ2tmanpDYOMxDGVpfnOJa4PKcxyfKcRy3LcxnD8lhZ2PNr2uhdZAGnKfvRu1z7M1unSX/WkU7ZtpoyO0hC6zlsIlaHVU8AtXQg04b0RrLtpDC3Yq3bJpUp07wBWqZC+wb4OBUQAJDqpwIGBqQ8lWUPM78FDCzH16mAAUBNUzEEUnJZprLs5FCmsjQvYSpt7Bg6P1U/dgR5ljFZN3jWsdow8MwW7tTFT5VVacMA0Lcv7fmpsk3yArB1O1Z+qmxcGwYAtq2NtvNTmwlXqADSKZJEgE68pEG+WrAcyRkgdiAqQKdiLCo/tU7Gmr25WKkipFd7joXnZs8anJ/E+1FEXExGFMSLrdCoLuokXjtRaEyvTOLDWIsTWzVkR+uHpOwCYBpVqwvFpTiJz6Mk7wobtE3iIRnZykrmVWxPdSE4VUCGxGVGVSYx2prYQzqJQMbEBXEKpDaqEKPcJoG28sxH0yQQV2aTgee9kjrGfT0L1NWLa82J95NIGSFw9EfqJFBX1l1uTqDobZTaaFgmUT9KhGZFJxn6Dh3r41YBdbRqvrjaJtEwarB2/CSQ1po4FhzbXKCtEpyRMqVwtq2uUjvOUNejI0C9jE3Ykk0mZUogif0zqYc2BTd/JmVOhP7YZ1KfOr2FrE/KpNRgJMBP6svI3iziCm+uo4bW0ZIyqW8jtSDWaVKmhFVA5UlFxiKuV55U1OjyXHtSJoX2qDypRGr2ypNKomavPCnzAUWXQu1JpXRsl0broXPWajN8l2Y57iBHHelJbUrAtrfLlPS+Sg2T2pRYZzXESTUeemstMylzb0FaoeLXo6pl6W+vXA8dpnJY9jC0zccpLNuYxcS5xU6eT/EO981s+V6nwH5e8g0+4gkM7rua8kmmwMaeixtYDxyEgdC2zm4kmaK/Zoz8FH3nCagOpIeaQNd8ClDkOOFEmaJnNDr1on/RNrPRfutPZDMDcxqnMEWWjtdD/6L4DpPvdYowS6DICWqw0tbcvSjBIB8NfzE8KEpvxA45Khu0kHo/RTFS3/slkHqqkM5BQwU/ReWc6V3SeUo6zzhF2KG5SxxpkU07d8m6zPqYu9TzoSRqXe75cB9zlxgCrUuXevG29MlKBw9IA1Y4iIF91qYYlkPPqgY79fp8TzHYwXuoaCfvoR7DYZM2JZbkvEZ4S2pymOoOX089GKV2Pfk9P/b+9D5OiSXZ2zGUU8ug3GdtSq0skL3L+t47OKVmXQfyOmXPFNipOGXfuXTqZE+fbdKm7OmxTdqUfefW4RhYhX2HGYoZ8tV32AxX6/I8bVPmTCG/T9uUjYTZCuxNsgr7CrQuZY4W22H9dRLsdXQRBDlZbBHOhft4LN2SY/Ehz1w7VI+afXfMneNUYTX24VBbCMzalFmSS/+oB8Wa+8cqyZCsvsNmeD745+nPRrN6/2yojUu3HTbD/XjrCyRrW7pH23DpvXcG9f05lwzXQlQOfVDgSDLcOXvMnvtw9GdbGjyWhRB2tKcyE8J52KfCqptRTDIVFt0MUjkd+muZ19NP98rM3hgKU5kZHFqywr13PbMuNNCoQ4UGzuOaZKq+UwKblKnORNCIRPWHOSdnRqfTnjqvxhnfqS6r0ShNZS3Oa8/aOAw2PGanyL2JeRFSmIIcpTNZoaRRwvmV14SQnJkc9e5WztGFIHudKudoX/ZWuB4ga+NAi2BybQHSPMB8OvTRmSpH6ELXeS1EktXWiXGFRs4glWdOeEZpXo3zKEx1Xo1LS32vAoESMuPcEMz2TLs6UizGuW9si4rIOA8R26KGPkh93Uz1mkCCcbDzApx62Zlqzwd1DYt80Atf02xryESUpWo8llFqNIrNe+y90Sh233xTjSYpzeskLgSbfsZrdMAuXqPTX3NAh47FAzq98IIOL0oLNgbMI9SP6KkiLs7zYh1Jh2Vk701gNFPn2qXGeUCnmjpG84BPNYGR5aIKSYbRNQhGM32sJkPO50HN/UDr0zbVLkf2SZuqCZIz4axdkOzYTjX3470vhanmo0OW9k2cvH5rH5+l1yZRHqp2kXLJK7MqZEGozKffglLpC3sZlNIXdqfYtfSFvXS0dGboMAoFyWaRsGvpC3sZCATNnsmKqn1h9+Ov1j5tS9HaZ60fGrXanPVTolbjgWZ0kDhtbfXJReacaaaNUe1bf2mjb/1DG7b15ybafPZYtZlSLyMwE+q5/Zk5sGOqzkR6BjpRXGZzptGH+Zxp9PKawwFkC3qm0faiNpPoGbim0IxbOyLRrOAGkZ7PYq9Tm8l0b6XNdHruSVvo9DzVbabTfTDbzDMsUOcY5pluM6GeM2ftmq0pkELcYOD7Z2oQS9Zx/0xNbQb7Gdo0jikgVkugqi3ziGhcpqZ5TN4k8IjQX8bYENQl8tI6xuoQPGOemrYxZpe8k1inFvwYIZ9OkkwtyBijM2Ti1IKOMTgjT2VqIYxRHeJ3RpMQx+gdGresUwtpDM2lhi5saiGPobrcVYQtlDEUh4hcaKGOITkE4EILDbVJzk5Knlr0YwguFyelTS3KGBRpHtm4RR2DiXxS09RiGLU5ZOFaphbjqNWUBA21RkJtVRLawqmhScuuFCeNFgoagtKcNFpAWHWVlSxTi21UdbXzdc30aM5EvzI1iGNzNaGTmhqksbqanUqYWkJ74Gp1KmlqKY6SXUN7VaeWTCBEEaYyNehick1niC1v6gnlnfB3qFw6ZGJzyx3K8FWulRliybiGnI+qBp0BugiVNrUcgATFhYGGjK1aA1PPDXNu7mBXXLRcOkgn69RypVvi6SUgClmHasLAgurToSpVKVMrYiDqMANR+ThUQh20ERLRGURL7NBeoAJvUMbqxHQNgGgXur7QwNJB2iUXYd4JWg5GG8oICBppatUbiKbAQOkgaJCL7tqUHh0MHQSNODUUcuTSfUAUsk6U7gPmDoIGYOkgaISpQSApDBqAzcAAGmFqrWOF7kl0aghRFe3jDPaxQt8kMjWoJNpp0ADsMxhAA7DPYAANdHC2oCSABqCtKNOuGmhLSiLru00NkT46QYvl0aWhtANmwVdglFEQjAVGq+4koo6w/DAKGujs1BfKo8FEI73AqcN1gVHLOIkNjTjli8EJdAzmVIHsOPWJfFOnSgIfYPFjm49iA8VAM2OQqyPK9K4rES/BQHCx3GggqARy04hKB0wMzIAZRNA4osBDn7iA1UDQsNw2os1jVFAxqre6ph0EFBqCnDXGVBUkMX80hsxYXVO3SmPI4HV1hukT3C46SKgTaJsSb4ZB1NjdbuZo4ILKCMyhX4ZMX/5GI2keAaxhtnFSaT/Ascww7SOEYVoyqwQwGsauU6zgE9DFzjD4BITDrhdFP+mDWQRMT9o1uQXCK6hnK6MRKmQZww0UXDzknfaBwQcCT/ud3ouHxNM+MO1D5GmfE4D2Y+Cw6HpMxiNGjpIOg09MHDQCPS+0HzPHkFTGnPZj4ZDqsLVfOcK6tYrxiI0DrsOMR/Icf2Z+y+CThMPR1KKZ9pNykHaY9mee2GxczEeKY6I+BzrtwxcXJxXBgfZTHjkjDbb2ypjBh3xrr46sHYPpf2ojJi70vqYqzx4NYodpL8uIzpr8SHvZDIYdpr0cRsYGC12kvRxHTIsG015OI20ZTHs5m6LbYPqXy2hae/gIa6+OGP+ATR2f29iw/gHTXvFYlTrM+HIseBCG3tEgBwMCjiXQovHMDJFRRBLiyDHSE2iT44FNU5uZCcVzQKBlJUGtVQgLzZJgzWIWoFkSrFnsduxRb8Y28RwUimmRBGaWo4Klawk0y2FhFkdKsLY4LliMWD7R83sODFYjJYRmOTJYftcJ6OSxWlKC8efYYIFZAohxcMQZD6wFnqODJWYlQIzDA+7MEhigRQdnCWBqSrh5PARMTR8yj4dRejtEGFMsFqBux8g8pkbLOUhs5VECTDlKIEuWAKYcJp3OYGokAYZgnlsSzAhklKA6gcyaGYitbQlYXRZLkCWYhQa7w/IOTCmcKWaIqA56ijXIVjMwthUzPdC16iCoZhBiIGrBUCpmEDIDenEQWCxC0lhvxUFgsQlJo1vFQWAxCwl2HGCzEqHGAaEF7poqg0HYzBB0CRYPGMkFfLJZbDANYbWkPAQY2xAmYYPBD+OQLeaKkVGwDmFkoT8QXAxEamu7ugp+6FKg+JhzaR/t3rzUIbhmIoLCN+8guFiJzOYMTPto96D4TVxhvNCn2EYQB8HFSqRQeGBrHwsja1odBFfQ70HhgRkflCq2TYKD4AoHDKatFlw3pOmotmuig+CKHTAs+OgguMIBY5sIph84jRpY/4gEwHnEhmkw7XPA2JbKDoIrHDCR9jKGWZHQfReYUwiuRD+qbbDiILiCSR4KDkz72ORtv1UHARbkCNtu1UGABUkimnTgEu0jS9juaw4CLEgTxqCYiCmCPGF6Bu8S44FEYbA6CLIgUxgcXKR9pAqDk4MgS0JGhbXIDoIsCS2pMWAOgiwJ7S2cSHMRfBI2WxhHMdOqJPADDg4CLZjsDY4OAi0J/MjPDgItCfyAC8ZskQR+wM0F8Mngh7TiHfRaMvgBq4NcSwY/4OgC+GTwA04OYi0Z/ICLC4xH7uPHGsJIKRn8vFNYY/DJ4AcsDsotGfyAg1PDB/yAMZtjagU/4Owg24JKxmCkIWDwIx9OHJj5bY41AxGX0ucXXxcBn5LmfPYgcJ5hOBJg5EXqw3UBs/6AMe7iSPT+qzN8yf7JHPJWH53v1vugK7c6N++xlNzqh9X4TdPomiL6ZdcU2lJd0+Za8K5h4g3qWgiuhehaSK6F7FqARlXXQnMteteiuBbVtRhci9G1iPCIjwi0rboWm2vJu5bEtaSupeBwN23YwFN2LUETq2upuZa9a1lcy+paDq5lHE2Sazm7lqGl1bXcXCvetSKuFXWtBNdKdK0gimbXSnGtQIeba9W7VsW1qq7V4FqNrlWE1OwatNTod3OtedewQTd1rQXXcHBpybWGCFtca9B9s6ujh4ESe9OEeCwQphn12GIgwR62wUN8PQ47HrLrobUeAuvhqz2k1UNPvUnFUFIP+fRwyR7C6eETPCTSwwl7iKOHInrYAw8t9BBAD9XzMLMeeuchch7e1UPevDkuQcg8HIGHhHnolof99FAsz5LzcJuexeahSp5l5rvzDjVYVR7FjIfm+EQN08rBVHrOeI/6xcNTek54DzfpOdo92hbPme7RtHgkEs+J7lGueI5yDyvpOcM9AoLn8PaoUDw8pOfo9nDqnjPbw4N7HAg8HKSHpfYc1R4G0nNGe1hHD+/pYRg9XKC3o92YBo5lD7PkOZA9fKPnKPaNGrCNnpPYwzB65hytBUwDXzg5MOfdq4k5RxHBec8XEj9zbowiWgYxDrH7QzDnxh8Kc94dA5hz85dCNcCpypc5UtAGc47cz/HIF20w58YQdo8s5tzYQWHOEdY5u/iiBnOOEM55hCKDNphzY/4QqDlY+KIGc24OW8bnIR5zLPBFDebclIHCnJt/nHm8CHNuDm6IsRBqvqjBnJuHnUmh3TOMOTf2zTzkjG8zodAYtu5Cwpwbu2bSmfFpJndJqX/405/+5P5W3pzp57w5+32En7k/Mnts4qZ9eUHBxbO5p9w7vODWRQjAo8sOt+85XN8AEe5dnF4B+XR9+f3m8uhKSU84euWccLgm8dXm9f7+xTPcrnFABeyZfvCMgeWfb59drNCKdPjo/WR/tMObPZuH8fr19tQ3/P5+f5903Lefbn/YXm13F1erMQktknP0wk/W322Way60Z3BvIdKCwZ+fnV1t7P4JdHZOPKAdDO/tk+8/2Vw84/aMHzxOzjYHS1XrC+7mt6vtzxd380ORpXkckx//anpo0/kX9PBffzU9PEzQnzmH768vj24Pvb++XBaFuWBzH4ttev7wy74nPrhc/9ivXXT485f76yseHZhveXRgvujx+cv9B93fvt8gwwmeTWTb6POX+5kk0InPX+4/shtXc9GPtvO9gBMPegpY4tPtnmtjC/zVbnduDvQk9CsqD3YX+92ry6v57sL9/YzOLYp5f79nExuR+hlaoH8mMWCv0P1+t4hOAnEjYjDow4unH15e7uYLW+xsA604TX306uLJTBbIBDyiYoDzFJLLlZe5MP0HnAuz7wGP5vuTzbPNxdPjezVg11OPKCwvuk5c2l5u2vGKhSwcCjJ+81J0q7OPuQuxubpFw+fUL1+un3AXwNo+XH876sPh7tucBo6HcjexORRdkpeit5q2crf7fZR4fano4+0VC/IYH5J434xO9nR7Kbc03EdnKbqkzgVvYUOpT7cX2xevXvzr5nJ3fdWDjBs3Eo3E91svDy83Z5vL//3JdemefjRwPeG4m2B6nHrdz576webs49WYPJN2SPl6NZabKY9W2AGOijyeEx6ujxffw/WNtUXjh6Trli3p9Cbmw/XTG11n7B6un57e7Hy4fnrH5c6H66cs9kfXwzOnPL6Rwtk431Wiwe2T7+ebSg/XL/sVykcz0TgkPF5hnlydffnkcrO5+Gj9xKgP6EHWjoYfkH1wtGxJOp6PpdbR/qEI4PXuIWVZQH1Z9TKXL6CiqzxEu0hDotGMVmdcjMQgnHbUvl7uO1HyYwCxfuwvty8/2DzZvlifXx0uFRlJnrkdPbAMR72zAre6Z2nH/WORWOJRBw/w4ajqTAoMSL/xe6h1IIFWp3enHUrTn+X9dKcTn93vLi42l1/QPUqy1ey1V6vxGy4S3bvHl+g9Qei/F+/VBc734j3JQOopEOci15AlHaerv1n8kFf9nFf9vejvxXvRSlrbf6OvP0BcNuunm0tOa7szZcN2gD7a7j9aFk2aF43dxGIWDzm2lmzUnqzPrTKz/3922wsSF0bgwfrlMfjV9sWBnSy1NqnRKMbvXqyfbXjRgcA/WF88Pd98/Xx79f3m8ov1xbP5jnJPf3/3ek7rs9dTDZOja5r/st2dby+W1Pn6Yi/6YHv55Pw2tZ+zuF8K0kcH4CNY7g9fv3x0zPYsiY+PEx/fVXJJvFGSgp+uX3+wfWZ36FmEn1/un+8erF9sLtcz9TkV1x6un/5i9++YsTfev3u4fvomeQ1JihE6oliAC/GZc4/G79bJeCqFnZLrO2j1bxfx35OhJM+flGpBl7iyi/nv5cHnFkVL0lpyZIK4qJ/uKE5yHLLGgjKlasgsvZvX+CXnAY2K16wlFNbJ0a1+SW2QRG6pivLyOnc1Nj+0GFqrueaAEuX2jX/1dchKiRICqkErsEQMUB+G5smtoUU7opbgAKvxvTCUoNy3TylJsdybgQHy0FJrrdUYY8D28nNxAqoOvjKYAQUTmN4MGqBBBrUCpqKz7AOemoeUqYzG/jprNb53O+s6JgAxBd6ToWpOwedckjRh8OegAHP2yQwv2V0S8K2lUkqpXrWcxhxIaYCwFvUxZBRqh9cTguCuhTK/vock0Fxqi5ol5FDVpvYQg2DwnrFtrZWSY9XEy6/DCgzpRrYtmjlMAVEvbNl6rxqS2KT/TF6PTeCHZqPvvY+xdCVFD1Tgh3iScyPkwZ2o3giCcAe2h6AIdzV8CJJwV9tz0IT3GK/TqAlM8Rw1IdYhpaStRp9sxf2ZIRRo4LcQCocINX92CIVbkvy71z8EXfj2VtiFDr9r4IVf7qL+zM78Xe7ps+reyCeYOPwmTgEK8ctyCsSjsX12HJCGHcbBdrK9joUMetHhWag8Za7+Ss7qo99/9oCR6PGOfnbUPpYP3jRmyQJnzAvujEA9377+NqWzxqGRpCEhdJZ9Na7u3/tit3t4vtvf253d+7fV638zhv5oyLva4a9kz/58JflnTzbniIkYfFZnb1ByzxruN4WUmofgtWlXj/r8GuBXrRYH/++2aPQxja3OHlmgHVvWjyy6jrFkvQSxYM4+2l5ezQrOT9bLE8G3dBFuX2w+2F69PF8fBa5B2jrIa4yZWQKuA8l8unv6yfq7Gf4ZTf27TdNPt6fJ4phc3fu3/3Xvf9zzg8/3/uc/1LwddAIH/cJsBfkZc8bRvB1mrc/fEh3p7z9pbzA+vNuk/fH2pP1DTRGDbLqJX/kU3bKXHBlTvG3/Dy/2l1s0Tyywr65evfhx3tU8Ljuc5yVylj3PGZ+uX6MRvuPwm3XFnJfdhPLZ7vLFoohk2c72jh5v6uzL3pjRmxvh+zAdEM3v1Ez57vxTDyl1raxcQkxZa+9vL9A1fnh5+TlBr0AN+PMfNpdn57sf54G5f3kJdbvpkpPx1TGXHDtpCNxkjjk3OYFf3uz9ko3y9nP+i90O5TyF7zrlzNvkENHr+beXm/X5dv/TB+v9miq29lfjilfsnl2uX3DIH5X59uX5bv/tt3YeHh35bJu/4sBn+H/+wIeru20V/+zlbntBMMj5RPu5ifLdd+o9P7Ti+Ba+a+K7NPsOfGfLzZabLDdabrTcYLnBctVyxXLFcr3lenFmlfeD8CTF+YHAJoM254eQnB8iGQQ+GbiINyQyMhkERhkKGZWMSgYeLrxUeL/QlNCqgICAi4CWgKGArIC30AWhNzLgnEEfhe4KPRcGQbG368Clm4F4OgM3WwbupAw4zA14vQ3c1xjwwRu4DTFwcWEgJMpAOKWBYDQDQXaGRkZrLmAHCvQ80PNAzwM9D/Q80PNAzwM9D/Q80PNAzwM9D/Q80PNAzyM9j/Q80vNIzyM9j/Q80vNIzyM9j/Q80vNIzyM9j/Q80vNIz+PQisW3e/wOi4ULSLgzOsLk4ENoDwS3ssTqcDwsjss2eDwSlgXXweZCjyJTrZZSkftIXMLpcU4czkS4uds/nG9ccTgBOlF7qk5cdMUFl5266KpTF1xwaumkAol9xP56y+GbgEs9rT9bbANL4+n673Wud2Jj8ksQ3QPdOMQOPDkSDqfHckbcoPwfovx5++R0bxH8XbxJuB8+er599vwvqIcy6e3VvgkDjmktxtBCxjkqDjnF5DW0FGNsmdUXU/W5tRY0BErUHLOo+lBjUluQ2RfJWXGVCy4NPmcvsamElGzZntQ4aeSkxK1W06CpFXz/SwsxJ5eGGHNSjSnFkH1zaZAcShENNaVUs0tD9Spao28hxaouDSGFFiTlWJNdax1aEZSA6EmDUiWFknIKmqpPpbk8SCiSAy54VbPc9Y6sRUIKUqRmT/eLpOAlVyk1WJVbqMfh5oCdjPnJCIYht5xQr9aqUZQNn3MIMeei5uV5xzTFEGMsrapKr3IrAUpUqsdtM/pQjSBVlRpqayEamWMcs4bgYyiQs9vwrcUTBs0plRa9b8nu8g2pSiw+puxzNiLcYogZF1BNOFGe9AximjSn1mLw+KW+/R06lFxDzT7lmLUoR4wWjVnVt5Lt3GgSQ2hSs6pyDfU2YjLkoMXHXDPel3ZclNaChNIUNTEnU5MqdEWqOV6evvQWHndUudXKSbMnVd7el5N36HCC+tuavaOVW1VOSvwlCaeVTgbxjjKnLd2RdHvK/5IEO0ZOZ+COdx0dOHdlnyCzlL+V0WPSPn43Gv9NHDTnFlvJuHNzl3womBdy9BF/WyhcTkmqtCJS8axP0IkW1OeCYYYqNScNMbeUK+7reZAcg0af2UBa76py0srJS281m4fA8ZE9lq/so8tDCtXHHIrmXLiZMmj2GkONOSWsL3motYr3KfhQk4jLA6YKaTlinanJlQFjVUxQdAh+hvLSRKilxtxcGVRqrCFVHwtu+yevyNDmkkL03hfoUeaMSDmXHJNww+QE85MBOxn0kxGMQ6UTvqSUaxVe6nNO2JKCKB7xd8xTimgkVVsyz+d0muArWsuUo5dsx6z6UKSFioWvwhPmVlIJjLFvxlueJNxaQHGIkpP6pEEyhxfMqtfqzegVuVEz4Gnuo4TCKMFv3uoc7Gmr0nwqoeEi/vZ3hKHFJFGStJrBVAdtDaqfGpZU3iFZC4eq0pl4+lIdoMY11aC1BgqkGnKUFFOrkVvNQw6lppY1xJQMr5N3nqBxq8rtRk7hkwpv68jtV9yB982O3a5wMlbvUOAWUne84lYJGSSWWmLMcIFcw79d4mSG3l7l5B1/QZWTd7w94aSVtyecrJ2TKreb/YNb/SDEur41bqsLWY3FrX7Q1fjN7Uonbzlt5+SFb084aeV2FU6db39Cr2lq7NV4t4L520PQ+m97oPqYbhJYt/r2Ynf54l/W54Tr9ga+2F6sz9/fXswuK2irSf9y++zF2jRK324WFV1X2X17uf7xoLVDb//t7nL7jLd8bTHKTVV2fHXsWE81qz/+W0m0v7it7fydFXAPXl3+8EYDZdIUju1tF30tfIuq9JvXfzjWwj283PETJ9vdBWo4co50btnUZm/WuuEN8bNeUG/XuqHW+xmtG/4rP6t1y7PWTQafuWdnD+4905+1aBq4at9qGrhs33xly0yWmSwzWqYlB/tWy1TLFMu01xJQbyDi3eAr6jdAIM1o33girPEQSUfllkjPpGfSiV098FBJb6Rzp4zXyoDmDNUNreKng+qNdPuQzoUrsBY6IINp30jnmtvAbTZ6zRlBIJyBEDYDwWcGNEqDKd1IJ4jKQOSTgQAlA0mZr0J6IZ2A7wNquIZg6JWvyknD8Q6kpKOGoreB3obB5FLSM+mZ9EI6D5X0Rnqr6NvQXiDnCnoNRYxW0k2rQV40fQfpiXTTh3QFB3wP6ZX0Fl3iL/3ly/v8du1bdtc69W+kDHDkLYbafLFIT7cT6hA15BJ89iVyXxFXsQpjX1vOqEXbUIuUpLG2zAU49QjoGqLPUjPDKoMKTJGEpna/V4YaQkoaY8mxECxpiMgFXhK/KaM6tNhariV7mEqnYYihpqiVy54WPHioJdWknTFDiTgoWoLSkPDRJQ6pRtEcizR6pnGosZXia1TBhq5p8CmkUnPydU6QFuBiS03cHE6D5takwfdzCZSE6ktNyiVL1AhpUM3JK91Ti6E+cDuzthQRdwytppFxkorHE2jk0qIiGqXMaKFmDTkVn4QPCd63nNVrrZXYUWhyvZQUW4UFQgFcg2+hxJC4ba1DaiVJVQlomEjwQVIqEmsomQBEQ9QYQ65wxsRi9kOpPqCUiCVYWKABZUhKPldfuOzZhtCazwg+0e7g1iEzeDP376QMjYVRTcFmoWlgbH0NORR67yQPodErqSWY6iINWUKQEIKWij4ELj7mWEtMicvJEgefavLVpxq4oBkGrTWnEIO0xP1aHZJWWOucsoVXkKFUYeYl4CPoRAaffCoRHZXapdsh4PmGd2Up3MoeMvqVwu3PmnN1dWigY6q1HD0JIaB8rOrh2qMrQ5GSkdFE0em4MogWRJbC3lBkq5RTKbEIo1gR+nxAZqiSSo2Cyg6RMYumwKK3fZtCaAG1lFaixQ3ZF421JB9CSEYhassIVwXfP6hI9V60hBBb4CIt+6I0lEmhVK71okkFlJRitQBRQ2nK9WrlAisBxIbkc/LRN7ax0UbWhoZAO4HlMngPRahFcdk0WsuOSDHn1Ozu9JAqI5VriMkusw7Rp1pbLKGi1ISARy1eJUcuMvMO31JWlLwsB7O9ZBQANbZQY+j6sph8EU1eakqZMgVNqTbNjbngSDEROTDUISXhCElZYo0aSpZCSAkc+krUlrK0Krli8YmqPkkrIbeihfeEHGmn5lSkREsRdr8mjY3lYcdaK/wcVFBpXNHG1iTSUkXPIUyiHyAv6ptE9aVYW5JrrDm3WqRyo5pzMmgpKdScpCbwEUW9GLTmKIRAYP9p5SYwt9EzP/gw+JrZ9DEUVj8npi/BSwks29xSwIrlsyhrRFk1bDvMYxI0No0VnbANkY+0lGPShiqW9nxIBX1MbhHVgzWoLecSMxe60bxTSmOqWVJCe83VY5J8UM05JImN0AhcZ8slJQ0+xyrJWykJpTCviNhcYKeUBC+stBJ8EUPVN40pFx+j1ERcBTw9C1qbGGuLVcqclgXtMD7GKREDg3Kx1ViyRGVZMareRx+ytMaleJTwlhZQfKu0GnJfVP7tSf8wJpubEtIJy3oiHZ0cUjelo18L2/5kd3m5ebLfPP3zGfeA08VvjPt7ds2Lb75+Y9w5LP/rMO6mr8OyVjFuEmYWxeNxQhqqNGxF2Tg8+LUIYa+BWFAEehpaCL6gbse2CJdYKgEuQlPJ2N1xJuCwkVB8QEiKQ2s1tQYPY4x8HkJsMRsXJNlpGZK01GqDOSMuYx2STxhiMdsSArINwUdshV5a8hj0hwbjmyIqaYIb+SHEjLGv1ZBaIAH1aW4Cn4grA44kgZMJ9h35zMPsR48d0JdmJUJonJ5etBZcCJAoQmvZuLhGQlQ03TXUNGNaSkvKGQ3zAepeU2zYaUvBI6MMUjJXVWDPEELyoDGnosaLE2wzDRobxuEgEV244lJRcspRq8eMDfPuU2ja0JqbaNRqROTQBKugfqixorWHobVYlUPxHF/YPqPFMxswmQi/vNcQIKRggoatTE0L8VyYlFSzz7UiBDhJA0pszYEbPkQniQOWlpIxW1uMmzBoLB4Wi2mCzeZKQ2y5cS+D2CVDaMg+yFA9xs2QYvWtZmEIxbWhcBwTRyZKIygRx2xNGugqkWTqEHxrCCVodGtxZYDbEG3GLIWAaQUTi5QYpRBsDXuNTyRqyphDjatWLcW3kAmrlIZUSsoltgTjOdvfa461SmipIU6XWrGcw2AR2SkOMVYvIWHRIW5OHCR6n2KDtbb4WAMsD0p8zAAZ55nsW46YeWCkzOGmSk6+eZ9yQs6XEln0mJdZxTrA/sVStGmKBCwdqqRmvJVnX6J2YMKa97kW2zxDkqTJJAA4WZQWGMI1EB6IOF2KjJti9B6pw3QVEgM+DrHl4FnWg4cGaAyqzSILDS1gS4fRhRtFUxJiyDmbSp3FMJQoISDDRkYS9UpSyfCEvpbu85SNjU8F6wvMPYIll4PgIImrM0SPHBJFMadQImjGBKKMibH/mszmlZtoYFEOUprHc6LEYuodIhLFkFKIqRD/aEBGgUL5HPGbwLksYR1vJmQkONYm6mPRmjzmIxRLiSBBMed+YQ0Jomlqgd0jhblEpmB5ZMUZINnBmGu35qWIJIqQESSXENhDkAiEjMru0Viqz8Tmwbss16App5Y84Yj8gFBmSogcW+X32IaYpOWA/Brg/kmRoBAyZsDUZ6HkWmOpQWFlKRKimv0p+yzsXz8g64WGTIj8TMe1Cn4BVbDDmS5OU2T9NE1VCRLoB9UWECZSM5aYFLQ4/MAty8zUeKglIgw3DiY2FoIBEMsqy5ffrMCWgKhVtDWpilQggfhOIFBKMtmJuEapYOLVnDMDRggu7sRlTRIa/fINEZB1J1qIJ4aI4wnlFEts0A1j/9+e9Ktlzk9Orl8nc/5yd/7T2/ny9zeXF1f7zfbiHuVvqdejhSA44tIXpn222PymXreNx57lY2SBp9/U678i9XocoumONWtCuwbrcTuhFU0tB9TjaJDTEJOvRbPUEKo511S0bhUerxGLLw+qKEIkBVReASZJ4fFD43RMJJRSG+q+xP+4qHA/tdRqauuID4vEVPCsSFWgxIV7wqhRggQOKFcGU/codB7+wZUheo0++Fo4kGHW0KmVFDUmj0BRONtS8ZyrNUQcYUJpKXJO4HdClaCaoq+izaeEGlUjnJtkKZXbySDW0PmHAptHAt4lHkdGzj861wgw6PE2SloZD5gsrLhSheikOOTgLYCgEqJU/G9MPYWDjRLzNg8B7hsFUU3cesa/MqLAbaJZYSIxZXBOVw5jvJHMhwPNdVXQSHCEUkIS7tsLCehQU0oeL6DCTNJgSLWlwNXg1D1nYFVQA+LTQmuweJojoYTjUGLmMAxctA9mBEIayDEX/B9t+cDVBY8CksCygxbTXWtpeJXAynYPGZOmYG29DyGr9bc7U1ZUp02JFJAzKmUkQPjYZrwPHuNBg8aIOSGy2xL+MAynwmajc87wJ5k2CZcaBtxW4RckCPFWMY5V3DfhuYgIEAYtGUkx+Ei4SjzXmUIT0oJtafU1WgwALdVc2LHFBGImRrUCMM0+lEjPzBdecY31JSY2Cx3BQ1Q9TFhjs9Fo8gGvIPTqJMDmtuYj2mQPGtEnbaZ/V5w8wxBbTDj2akS5b4y/T9h2KkpoSqAUTzVX9MjZLHzYJwT3H61mzPOIxzEQfxNrVhzgUTOqUpyGmbXoo2IfUhzomCWiN+QcU8SL2LybcvEFnjFp5idgBnYlM+9xJIosNyWwKLp/T6zaNJQQEbdC1kgg0IwwhTMwliGiuOKGRnyC1HKSZJSkoqOtKLBDSea5hu+yOQ5j8CoDUhimBK8ZvXUdVDRE4kw0ZBlXhxKqNpTRrEqHdI+HX0Faq4SoHZpXgko0TGQW6hSf6pghbqiPES01iCaYVZ8tYubQIJi1tlDMkKSIKAUblzaL6hmGVAShusDDcwFgiCnWoD4xsvwWxRCIetBKkWbxYTFGoX0WrDBYpwpVmEKPqzVB2c0ZN6SEz7lFuoaOIutXwY6GgsTHUpCmctbKb8PQSipFgkeg1TBU1VjBhYAAKARCSCUHgWhgfs5DCxUXR5h0pMA65MBulZa0dm1HhMokdCbE0Q3cGIkJ827BvuEClmKcB5MoclHIQyxQEhwo8aUMFXNVK94jKhVxEeMLLsQh1GqGJUyKJUeJGbtOdhFVRmADY4vBCfttCb8KHvyEBT1Rl5/082/Mkd9/9Lsvv7z/6Ydvdnb5p8vN2WqU9otr6r87/+yDd73N/nD9g8X4As277rn9ggGD8JMhBgA6hCowOkEt8G2PhJNwEQk+h5YDuiRK62rE5i3VS0oZzsECl/RAOJDV6KtWT7hu3k23ifKhhDYWmBwzc/L+OYurAvwJLcIkWRtzFiyUaWh81hwQSx4vWa310D0Rf/sCZsfRCejX7y6IYgkWXz5fP939eBRf6cHuknBl66fbV1f9huThos5houb78ntuMDMTxOu6O1Dtktmj1Cq43A7RdohSex1Ddo4W9snu4tmGO+gmk326vny2vbDQNFxaR0T78+6J3g7Jsd5v7ELjXesIxOdeWgetp28KgvCu3V9G4zho3aH3ZM6B/+aIf/TwaGKOZNLl3uUsk3YR9cV236/c9ruvrIjPzVXwcKcZbhN21JTeKf7pD2614xrsNxevzs+J2PyHk5jNvwUBO6w6Bv0kAOSjM4sNyww8vn683x9tl19ZaFniElHmBnS/Q1S22Ca33faIyHFrERzKLSHQevzFQ9EleXkl1Y/Do/HKY2JwKHgjXolRJ4sgpZ2SIH6pBZkkGAqXZ3EG8UgMJeBm0okfLNNRhU6UVmMZ/HFxkm9GAStcTqtG5ohHPxdYImDlChOOCR3vI4tEtcQIswhjSDOtQZtrpT+P97slEtb3qzFGuNh3jAJ2HCOMAFt1gGHEZaxkftLJgnwdvfwXjQIWc/mbRAFDF5oIb4+szwAtYb56kLDjGbOZXLJ7FLAakDwx0QQsWCviUy2BvL5fjYJiIOIRV6Ri+Di83qKA3bVSltq2kDzKffjJooqFi7XxznHA7goT9m6RwOJJlLC/RySwu/D9q2OBsVtPQ4FZONw5Flge2FF/ZgAwFvtvAcD+UwOAEQ/2EJrhiXT4twBgd3P+UI5bh+WbOSY79m4G0T4NFQrzz+Y6jv9l2wqGvu+p20fpjVP0FxeUWA5vjwTy1cfyESXvYmsZpJmtvVq/eEkU3A+2V/sPX+8vb8QBof7RYL7NdemdhvPnI33wipM7B3NoLxMbftnIXnTwv0hgL7tU0uN6GePWI0TZ4688rtc/1CQxoP8VQkT9FsXLZvEo/vsSIA/W08LjHW0j0v7uu+jdo3jZouxhvK6f0UIdMl6vRhzGuE8guXlfC306xPXK3AYiAnKuyAicDyeBWk5Cufxdo3vd0tqgACRkl/1i0V1H3C+oAZyFa/gA3AX5owlfwGJhQIWw14RLIVRHU19NxEcJyF0P+8N9Ek0Wmtp4hlSsdG5mXGNxLZq+0sMML7oZYztmTV9Z4jS37GMm69FB0zfXah4jS29nUR2GGeNQI2aom0pAFsFfqAS8pIEjPmKZjHfVgolJFjd/huJaDfYzSsC+OP/zlYB3c05vVIEOPmRQtxJdD2q783oEFrDTJKC/sSIw5GLRtaXYBacTrd9fohqdLynf+2d+noI/rJE3j8qS+Y8yIH0j4TtZQ/5FBuTTzfri3r1/vncPA4L90Od/59H4cv/03gebH+7NIyK+C7dv3DX/2OvjhubcoT6//jIDpP1q2Bxh8vOXN39f7K6fPeFIO0TQXX5A5+H6csPP9/3HbHuL9U83W7am/lYBL49+g1aP7+bGOTpdjzenf4NferxavzCz0t9E1jWyNJtzeD4Wcy2E9F+pNfhNzP0lf9YRPqPz37+JuX+X3940wsv2OFgNfwtW3X+76r+IuuhXLubaj+ua+IqchbDXZdnVmMrgpbaSuO3dZbGDmKsaiB+BAazikPiPJ+YeOIr/lDO9R4olSAlxZNXizh5Hil0C9En3EuCHAr58crl9SfghO0U/3j57fr599nz/YHdxsXmyv/5Nw4+2rzdPOzE5W59fbf70/wFVdpxE').then(json => {\n",
       "   const obj = Core.parse(json);\n",
       "   Core.draw('root_plot_1779222972440', 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_1779222972440();\n",
       "</script>\n"
      ],
      "text/plain": [
       "<IPython.core.display.HTML object>"
      ]
     },
     "metadata": {},
     "output_type": "display_data"
    }
   ],
   "source": [
    "%jsroot on\n",
    "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
}
