{ "metadata": { "name": "" }, "nbformat": 3, "nbformat_minor": 0, "worksheets": [ { "cells": [ { "cell_type": "heading", "level": 1, "metadata": {}, "source": [ "Calculation of control fields for Hadamard gate on single qubit using L-BFGS-B algorithm" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Alexander Pitchford (agp1@aber.ac.uk)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Example to demonstrate using the control library to determine control\n", "pulses using the ctrlpulseoptim.optimize_pulse_unitary function.\n", "The (default) L-BFGS-B algorithm is used to optimise the pulse to\n", "minimise the fidelity error, which is equivalent maximising the fidelity\n", "to an optimal value of 1.\n", "\n", "The system in this example is a single qubit in a constant field in z\n", "with a variable control field in x\n", "The target evolution is the Hadamard gate irrespective of global phase\n", "\n", "The user can experiment with the timeslicing, by means of changing the\n", "number of timeslots and/or total time for the evolution.\n", "Different initial (starting) pulse types can be tried.\n", "The initial and final pulses are displayed in a plot\n", "\n", "An in depth discussion of using methods of this type can be found in [1]" ] }, { "cell_type": "code", "collapsed": false, "input": [ "%matplotlib inline\n", "import numpy as np\n", "import matplotlib.pyplot as plt\n", "import datetime" ], "language": "python", "metadata": {}, "outputs": [], "prompt_number": 1 }, { "cell_type": "code", "collapsed": false, "input": [ "from qutip import Qobj, identity, sigmax, sigmaz\n", "from qutip.qip import hadamard_transform\n", "import qutip.logging as logging\n", "logger = logging.get_logger()\n", "#Set this to None or logging.WARN for 'quiet' execution\n", "log_level = logging.INFO\n", "#QuTiP control modules\n", "import qutip.control.pulseoptim as cpo\n", "\n", "example_name = 'Hadamard'\n" ], "language": "python", "metadata": {}, "outputs": [], "prompt_number": 2 }, { "cell_type": "heading", "level": 3, "metadata": {}, "source": [ "Defining the physics" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The dynamics of the system are governed by the combined Hamiltonian:\n", "H(t) = H_d + sum(u1(t)*Hc1 + u2(t)*Hc2 + ....)\n", "That is the time-dependent Hamiltonian has a constant part (called here the drift) and time vary parts, which are the control Hamiltonians scaled by some functions u_j(t) known as control amplitudes\n", "In this case the drift is simply a rotation about z and the (time-varying) control is a rotation about x\n", "In theory this system is fully controllable (irrespective of global phase) and so any unitary target could be chosen; we have chosen the Hadamard gate." ] }, { "cell_type": "code", "collapsed": false, "input": [ "# Drift Hamiltonian\n", "H_d = sigmaz()\n", "# The (single) control Hamiltonian\n", "H_c = [sigmax()]\n", "# start point for the gate evolution\n", "U_0 = identity(2)\n", "# Target for the gate evolution Hadamard gate\n", "U_targ = hadamard_transform(1)" ], "language": "python", "metadata": {}, "outputs": [], "prompt_number": 3 }, { "cell_type": "heading", "level": 3, "metadata": {}, "source": [ "Defining the time evolution parameters" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "To solve the evolution the control amplitudes are considered constant within piecewise timeslots, hence the evolution during the timeslot can be calculated using U(t_k) = expm(-i*H(t_k)*dt). Combining these for all the timeslots gives the approximation to the evolution from the identity at t=0 to U(T) at the t=evo_time\n", "The number of timeslots and evo_time have to be chosen such that the timeslot durations (dt) are small compared with the dynamics of the system." ] }, { "cell_type": "code", "collapsed": false, "input": [ "# Number of time slots\n", "n_ts = 1000\n", "# Time allowed for the evolution\n", "evo_time = 10" ], "language": "python", "metadata": {}, "outputs": [], "prompt_number": 4 }, { "cell_type": "heading", "level": 3, "metadata": {}, "source": [ "Set the conditions which will cause the pulse optimisation to terminate" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "At each iteration the fidelity of the evolution is tested by comparaing the calculated evolution U(T) with the target U_targ. For unitary systems such as this one this is typically:\n", "f = normalise(overlap(U(T), U_targ))\n", "For details of the normalisation see [1] or the source code.\n", "The maximum fidelity (for a unitary system) calculated this way would be 1, and hence the error is calculated as fid_err = 1 - fidelity. As such the optimisation is considered completed when the fid_err falls below such a target value.\n", "\n", "In some cases the optimisation either gets stuck in some local minima, or the fid_err_targ is just not achievable, therefore some limits are set to the time/effort allowed to find a solution.\n", "\n", "The algorithm uses gradients to direct its search for the minimum fidelity error. If the sum of all the gradients falls below the min_grad, then it is assumed some local minima has been found." ] }, { "cell_type": "code", "collapsed": false, "input": [ "# Fidelity error target\n", "fid_err_targ = 1e-10\n", "# Maximum iterations for the optisation algorithm\n", "max_iter = 200\n", "# Maximum (elapsed) time allowed in seconds\n", "max_wall_time = 120\n", "# Minimum gradient (sum of gradients squared)\n", "# as this tends to 0 -> local minima has been found\n", "min_grad = 1e-20" ], "language": "python", "metadata": {}, "outputs": [], "prompt_number": 5 }, { "cell_type": "heading", "level": 3, "metadata": {}, "source": [ "Set the initial pulse type" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The control amplitudes must be set to some initial values. Typically these are just random values for each control in each timeslot. These do however result in erratic optimised pulses. For this example, a solution will be found for any initial pulse, and so it can be interesting to look at the other initial pulse alternatives." ] }, { "cell_type": "code", "collapsed": false, "input": [ "# pulse type alternatives: RND|ZERO|LIN|SINE|SQUARE|SAW|TRIANGLE|\n", "p_type = 'RND'" ], "language": "python", "metadata": {}, "outputs": [], "prompt_number": 6 }, { "cell_type": "heading", "level": 3, "metadata": {}, "source": [ "Give an extension for output files" ] }, { "cell_type": "code", "collapsed": false, "input": [ "#Set to None to suppress output files\n", "f_ext = \"{}_n_ts{}_ptype{}.txt\".format(example_name, n_ts, p_type)" ], "language": "python", "metadata": {}, "outputs": [], "prompt_number": 7 }, { "cell_type": "heading", "level": 3, "metadata": {}, "source": [ "Run the optimisation" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "In this step the L-BFGS-B algorithm is invoked. At each iteration the gradient of the fidelity error w.r.t. each control amplitude in each timeslot is calculated using an exact gradient method (see [1]). Using the gradients the algorithm will determine a set of piecewise control amplitudes that reduce the fidelity error. With repeated iterations an approximation of the Hessian matrix (the 2nd order differentials) is calculated, which enables a quasi 2nd order Newton method for finding a minima. The algorithm continues until one of the termination conditions defined above has been reached." ] }, { "cell_type": "code", "collapsed": false, "input": [ "result = cpo.optimize_pulse_unitary(H_d, H_c, U_0, U_targ, n_ts, evo_time, \n", " fid_err_targ=fid_err_targ, min_grad=min_grad, \n", " max_iter=max_iter, max_wall_time=max_wall_time, \n", " out_file_ext=f_ext, init_pulse_type=p_type, \n", " log_level=log_level, gen_stats=True)" ], "language": "python", "metadata": {}, "outputs": [ { "output_type": "stream", "stream": "stderr", "text": [ "INFO:qutip.control.pulseoptim:System configuration:\n", "Drift Hamiltonian:\n", "[[ 1.+0.j 0.+0.j]\n", " [ 0.+0.j -1.+0.j]]\n", "Control 1 Hamiltonian:\n", "[[ 0.+0.j 1.+0.j]\n", " [ 1.+0.j 0.+0.j]]\n", "Initial operator:\n", "[[ 1.+0.j 0.+0.j]\n", " [ 0.+0.j 1.+0.j]]\n", "Target operator:\n", "[[ 0.70710678+0.j 0.70710678+0.j]\n", " [ 0.70710678+0.j -0.70710678+0.j]]\n", "INFO:qutip.control.pulseoptim:Initial amplitudes output to file: ctrl_amps_initial_Hadamard_n_ts1000_ptypeRND.txt\n", "INFO:qutip.control.optimizer:Optimising pulse using L-BFGS-B\n", "INFO:qutip.control.pulseoptim:Final amplitudes output to file: ctrl_amps_final_Hadamard_n_ts1000_ptypeRND.txt\n" ] } ], "prompt_number": 8 }, { "cell_type": "heading", "level": 3, "metadata": {}, "source": [ "Report the results" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Firstly the performace statistics are reported, which gives a breadown of the processing times. The times given are those that are associated with calculating the fidelity and the gradients. Any remaining processing time can be assumed to be used by the optimisation algorithm (L-BFGS-B) itself. In this example it can be seen that the majority of time is spent calculating the propagators, i.e. exponentiating the combined Hamiltonian.\n", "\n", "The optimised U(T) is reported as the 'final evolution', which is essentially the string representation of the Qobj that holds the full time evolution at the point when the optimisation is terminated.\n", "\n", "The key information is in the summary (given) last. Here the final fidelity is reported and the reasonn for termination of the algorithm." ] }, { "cell_type": "code", "collapsed": false, "input": [ "result.stats.report()\n", "print(\"Final evolution\\n{}\\n\".format(result.evo_full_final))\n", "print(\"********* Summary *****************\")\n", "print(\"Final fidelity error {}\".format(result.fid_err))\n", "print(\"Final gradient normal {}\".format(result.grad_norm_final))\n", "print(\"Terminated due to {}\".format(result.termination_reason))\n", "print(\"Number of iterations {}\".format(result.num_iter))\n", "print(\"Completed in {} HH:MM:SS.US\".format(\n", " datetime.timedelta(seconds=result.wall_time)))" ], "language": "python", "metadata": {}, "outputs": [ { "output_type": "stream", "stream": "stdout", "text": [ "\n", "------------------------------------\n", "---- Control optimisation stats ----\n", "**** Timings (HH:MM:SS.US) ****\n", "Total wall time elapsed during optimisation: 0:00:03.691550\n", "Wall time computing Hamiltonians: 0:00:00.205678 (5.57%)\n", "Wall time computing propagators: 0:00:03.127814 (84.73%)\n", "Wall time computing forward propagation: 0:00:00.047317 (1.28%)\n", "Wall time computing onward propagation: 0:00:00.046361 (1.26%)\n", "Wall time computing gradient: 0:00:00.253621 (6.87%)\n", "\n", "**** Iterations and function calls ****\n", "Number of iterations: 8\n", "Number of fidelity function calls: 17\n", "Number of times fidelity is computed: 17\n", "Number of gradient function calls: 16\n", "Number of times gradients are computed: 16\n", "Number of times timeslot evolution is recomputed: 17\n", "\n", "**** Control amplitudes ****\n", "Number of control amplitude updates: 16\n", "Mean number of updates per iteration: 2.0\n", "Number of timeslot values changed: 16000\n", "Mean number of timeslot changes per update: 1000.0\n", "Number of amplitude values changed: 16000\n", "Mean number of amplitude changes per update: 1000.0\n", "------------------------------------\n", "Final evolution\n", "Quantum object: dims = [[2], [2]], shape = [2, 2], type = oper, isherm = False\n", "Qobj data =\n", "[[ -9.07307932e-06-0.70710599j -9.73151123e-06-0.70710757j]\n", " [ 9.73151123e-06-0.70710757j -9.07307932e-06+0.70710599j]]\n", "\n", "********* Summary *****************\n", "Final fidelity error 8.914791127523358e-11\n", "Final gradient normal 0.006605857936333281\n", "Terminated due to Goal achieved\n", "Number of iterations 8\n", "Completed in 0:00:03.691550 HH:MM:SS.US\n" ] } ], "prompt_number": 9 }, { "cell_type": "heading", "level": 3, "metadata": {}, "source": [ "Plot the initial and final amplitudes" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Here the (random) starting pulse is plotted along with the pulse (control amplitudes) that was found to produce the target gate evolution to within the specified error." ] }, { "cell_type": "code", "collapsed": false, "input": [ "t = result.time[:n_ts]\n", "\n", "fig1 = plt.figure()\n", "ax1 = fig1.add_subplot(2, 1, 1)\n", "ax1.set_title(\"Initial Control amps\")\n", "ax1.set_xlabel(\"Time\")\n", "ax1.set_ylabel(\"Control amplitude\")\n", "ax1.plot(t, result.initial_amps[:, 0])\n", "\n", "ax2 = fig1.add_subplot(2, 1, 2)\n", "ax2.set_title(\"Optimised Control Amplitudes\")\n", "ax2.set_xlabel(\"Time\")\n", "ax2.set_ylabel(\"Control amplitude\")\n", "ax2.plot(t, result.final_amps[:, 0])\n", "\n", "plt.show()" ], "language": "python", "metadata": {}, "outputs": [ { "metadata": {}, "output_type": "display_data", "png": "iVBORw0KGgoAAAANSUhEUgAAAYwAAAEZCAYAAACEkhK6AAAABHNCSVQICAgIfAhkiAAAAAlwSFlz\nAAALEgAACxIB0t1+/AAAIABJREFUeJzsfXe4VcXV/rvuvYAIAoooVdDYe4tGzRdRY2yx5LP/NMaO\n3VgSEz+Va2yJiV1jNMHeo7GjBk3AGruIhVhQQUBAEBCpl7t+f8we99pzpu1z9rkXZL/Pc56zy+yZ\n2bNn5p1VZoaYGSVKlChRokQIDe2dgRIlSpQosXSgJIwSJUqUKBGFkjBKlChRokQUSsIoUaJEiRJR\nKAmjRIkSJUpEoSSMEiVKlCgRhZIwSizxIKLhRPRzz/3rieicyLhGEtFRxeWu7UBEzUR0e3vno8Sy\ni5IwSrQLiOhTItopJiwz787MtyfPHU5Ezxn3j2fmCyOT5uTnytfaRPR3IppGRDOJaDQRnUZENbUV\nIhpMRBNqiQOefJco0RYoCaNEe8HbcbcHiOh7AF4G8BmADZm5B4D9AWwBYIU2SL8xFKTeeShRwoeS\nMEq0OxKp4Xki+iMRzSCicUS0q7g/koiOIqJ1AfwFwDZE9DURzUju30JEFyTHKxLRY0Q0NYnrUSLq\nF5mV8wE8z8xnMvMUAGDmD5j5UGaelcS/FxG9S0RfEdG/kzzpfH5KRGckUslMIrqHiDoRURcATwDo\nm+R7NhH1SVRM9xPR7UQ0C8AviKgvET1CRNOJ6EMiOjqyDHv43jspwwuI6IUkD48Q0cpEdCcRzSKi\nV4hooAjfSkQnE9HHibR1KRFRcm9NIhqVvOM0IronsnxLLOUoCaPEkoKtAIwF0BPApQCGiXsMgJl5\nLIAhAF5i5hWYeSV5Pzmm5NnVkt88ANdG5mEnAPe7bhLR2gDuAnAKgJUBDAfwKBE1iXzsD2AXAKsD\n2BjA4cz8DYBdAUxK8t2NmScnz+wF4O/M3D2J+x4A4wH0AbAfgIuJaIeIvDdEvPeBAA4F0A/A9wC8\nlDyzEoD3AQw1wu8DJV1tDmBvAEcm1y8A8GQigfUDcHVE/kp8B1ASRoklBZ8x8zBWi5vdBqAPEa1i\nCedSyxAAMPMMZn6Qmecz8xwAFwPYPjIPPQFM9tw/EMBjzPwMMy8G8CcAnQFsK8JczcxfMPNXAB4F\nsGkg3y8y8yPJca8krrOYeSEzjwbwNwCHhTIe8d4M4GZm/oSZZ0NJPB8w87+Sd/k7gM2MaP/AzDOZ\neQKAKwEcnFxfCGAQEfVL8vliKH8lvhsoCaPEkoIv9AEzz00Ou+aNhIiWJ6IbEvXQLACjAHTX6pQA\npgPo67nfB2r0r/PJACZAjbI1vhDH8xB+h8/FcV8AMxKJRGO8Eb8Vke89RRzPBzDVODfzKo3045GW\nza+hCPAVInqHiI4I5a/EdwMlYZRY2uAylOvrZwBYG8BWiZpne6jOLYYwngawr+f+JABSz08ABgCY\nGBG3Ld+m4X8SgJWISHbcqyFLKi7kfe8Yh4PVjOOJAMDMU5j5WGbuB6Ui/DMRrRERX4mlHCVhlFja\nMAVAfyLqIK7JjrEr1Mh+FhGthEq9PODuRIcC2DYx8K4KfGvgvZ2IugG4D8AeRLRjkv4ZUCPzGJXM\nFAA9k3is+UhUPy8CuCQxlm8MZTe4IyL+vO8dQ6BnJsb0AVB2m3sBgIj2J6L+SZiZUOTTGhFfiaUc\nJWGUWBJgc7F1jYCfAfAugC+IaKoIq8NfCWVX+BKq830iNm5mHgdgGwCDALxLRDOhjOCvApjDzB9A\nGY2vATANwB4A9mTmltB7JQb7uwGMS7yY+jje++Ak/UkA/gHgPGb+l+U9TeR975gyfxjA6wDeBPAY\nUkeELQH8h4i+TsKcwsyfOvJV4jsEas8NlIjoJqhGN5WZN3KEuRrAbgDmQnmcvNmGWSxRYpkEEbUC\nWDMh0RIlALS/hHEzlLuhFUS0O1SlXQvAsQCub6uMlShRokSJLNqVMJj5OQBfeYLsBeDWJOzLAHpo\n3XKJEiXqiiVqFn6JJQNN4SDtin7IuvZ9DqA/su6BJUqUKBjMHFqmpMQyiPZWScXA9OYoRz4lSpQo\n0Q5Y0iWMiVB+7hr9YfF5J6KSREqUKFGiCjBz9KKWS7qE8QiSZRGI6AcAZupF4Uwwc1W/b75hEKXn\nAOPwwxlz5qjjo46qfGb+fMbEie44AcYXX1SXH1d8zzwTDrfGGowuXYaCmTFhAmPNNatLa7313PE/\n9JA63nJLxsMPq+Mvv2S8/XZx7+vL21tvxYcfOnRo5tnHH7fHee+96njOHMbs2bXl79Zbi3vXvfeu\n7fn337eXRZ7f6NEqLhmvrgO23wMP1F6GvXvXp/60tDAWL66+LObOZSxY4M73E08Uk8833kjL3Cx/\n83fOOe7y+uYb/7Mq3/nQroRBRHdD+YyvQ0QTiOhIIhpCREMAgJmHQ/mtfwTgBgAnFJ2HiRMBZvV7\n7DF17ZZbgN/+Vh03NgJHHw08/XT6zG9+A/TrBzz8MDBmjD3eKr6FF8zAtGn+MJ99BnyTLCrx9tvA\nRx/5w0+YANgWzHAtojFuHDBihDp+7TX1/gBw5JHAxhv707Jh0iTggQfyPdNaw/SwFsdsiQ7JFMDB\ng4ENN0yvT5wIdOyYL41XX3WXX148/DCweHFc2PffB6bUwbI3f776X7QovTZnjjv8vvsCt91WW5pF\ntx2Nnj2BX/6y+udXXx049FD3/dhvFYKMJ1QWb78NfPFF5fX77wd+7txyrHq0t5fUwczcl5k7MvMA\nZr6JmW9g5htEmJOYeU1m3oSZ3yg6D7Nnq/8ZM4A990yvf/21+m9sBIYNA26+Ob03aZL632cfYOut\n7fHWo9Kvsgrw/PPu+zLNmPQnRixo8dln2fOFC9Nj3ZnMnYuqcMklwH775XvG917jxwMnVDGk0ITx\n3nsqDo2PPsp2lDH48sv86fvgIjkT668P/OIXxaYNpIRp1gMfmpZQRfesWcDrr1def+ghNQBy4b//\nBfbfXxHyrFnucA0F9aayjusB0l13Ab/+dTYcEfDEE/Y4hg0D/vGPYvIjsaSrpOoOXQHMjkif2yq/\nbsQdOgDz5tnjrYUwzjsPuPRS+z3ZodkxGEDcSDwUZuJEYNCg7DXZgep3r/Zdu+ZeWlDleeFC+2ju\n0UeB68VMncGDB0fFqTtFXR6aFGMkhYkT840sL75YjVTrgeWWy57L/MeWBQDceitwT7LDhS4TPTiI\nQV6pzESe+jR6tBrs5YFZFj/7GXCEZ/nEhx9WI3bA/+189YUZeCMw3D3wQGDs2Gy71Md/+APwxz9W\nPpN3QFMrlnnCmDlT/bsIo9FwLnzuuZS5fRW7FtXJ9dcDZ51lv/eVZdbKU0+pTk7lZ3Awbxqujk4/\na+skZAV98EH1H3rXgw/OSiYaJmGceCJw992V4VpbU7GbGVhjDeCQQ/xpAvkJQ793p07ukZuJ/v2z\nJBXCv/8NfPppfPg8nWevXtnzagnj8MPVD0jrSEx9/jxZIrFDB3+4IrHppsBxx+V7Jk9ZANm6a/YH\nEj7CePppYIst/Oncd59Si9skjBDmznWrx4vEMksYLS3AueeGJQxdQXRluPbayjA21CJhbGRdJEXh\npJMqr+26qyKxvCqpakjNNqIJpXXPPcAzz1ReNwnjz3/Olq/GtdcCffqo48WL1ag+NFrLA93ByfL4\n5JNsGN87Tp9eTD7mzKlND96zZ/a8FluKfl9dJqecon4+7Lij+vcRxmefxavZYlFU+btgG+zYIFVS\nCxakUkmeOLQ9VcPVRv/v/7Ln77zjtiMWqR4PEgYRrUNEzxDRu8n5xkR0TnFZaB9MngxceGE6is6j\nkjLD2ODrjF0qFY2VV3bfc8FUjRWhkrLBVvFjiFNLchI2ldSLL1aWjVTD6TxX0xm6vrGNMOTxxx/H\n66dtZdHS4jeWaqywQmVHkAdduyr1zKOP5n92xgy7fl6Xw6hRwDXXqGPX99bfzUcYgwYBN9zgvi8R\nSwQ+I3wIzc3q31efHnnEfU9CxvH448rukRfM7noo8fjj2XOf1KpJ6JBDaiePmGbwVwBnQ+2yBQBj\nkO68tdTC7ABiVVK+ZyQeeADo69iKZ8AA5Vnkgh4p5unQpfqIqDaVlA8LFwI7GBuGxhCGLUyXLvZ7\nptfHggWV8dk68LyNQY90dUM3G6q+nseQbcvDV18Bd94Z9/ybxtKaed6poQG47DJgr73in9FYay3l\nJWamK+tIiKT19wwZvbWjiQ8zZ8YPnPISxqhRaR0+//xw+NGj02NfGcg6Wa1019pqlzB0fDNnqgGU\nOXDTTjo2aBK6665sW6oGMYSxPKt1nJLEmQG0samleJgftGjC+Pe/lRRjw9SpwFtvuZ/VDU6O+HR+\nV1zR/oxpbyhCJSXLSHvJLFoEjBwZn5aPMPR7hiqxvF+LhGFCNzpbHhcvTtMoQpqJhdnw88Rj5jNP\nvmfMyHpCmSqpmPi0ZBEKJ43iCxdWuroy5yOBPGFffBG44AJVh7XEBMSXlS+cvGcOaGLjD0kYBxwA\nDBxYSRhmOFM9rYnfRywxiCGMaUS0pj4hov3g3/d4qUAsYeR1EbSNzGzwqTh0HHIkFuo4qiGMPBKG\n9um25dtHPDFkEsqHTcKIbYBvvlnptTN+vPrpeG15lO+Ux11SPyfjzNNx10oYMi0z3aOPVt42eZCn\nHPT90PeUKquJE4GrrqoMo1WsMVJ2XglD29NCNpm88EkYsd8xZMOYmuwAYw6yQm1Q369FfQfEEcZJ\nUJPm1iWiSQBOA3B8bcm2P/QHdY2AXRJGqPHrD1MEYeTpLEwbhnx23jzgyisrn8mj8tIumzaJS6dl\na/i+d3HdM8u4FgnjzTcrDfXrr68cC/Qo7d13K6UJ+f2uuy6cjp5UqZ8LeaC5UAthhDBsWLxHV0jC\nOPNMYMiQ7DOxhBHjdqs7thi30Vo7QaAYiZVITW5lrj6+kIShv0tIwjCh79ddwmDmj5l5JwArA1iH\nmbfj78DuWkWopGzQz9k8Qe69N023aMLwSRgvvgicdlrlM3kkjOWXV/8+CcPUv8t85CEME/LdfDYM\nWyO1eY59843qZHSjGye2CNLxygZ4yy3+/AHpfJU8bqg2+FQLIYTscrb4X3mlcv6GfNZlw7j2WuDG\nG7PP6LayeLEaoJx6avb+979vz6ct7RBhPP88cM45lXmsBWPGAAcdVJmXWBAB3/ueWg2h2kl8JmGY\n6ZtzhQCVVkjCKEol5VS4ENEZMk1xPckEX15b0u0LU8JwNVStknKNGFzE869/VYaVrqBFEcaHH6p/\nkzBiOqw8Da1zZ/VfLYH6Oq/Qe0oJQ+c5j05YQndARHbbiTSA+9KYNw/41a/Scz3b3SZh5Bltmh5G\ntdgwYgjj9dft5WBrF7LO2urOttsCL72k7l12mZqXsf/+wA9/qO7r2dSys3OVTYgwrrgibj5UHjz0\nkBrUnXCCcjzZYIN04Lfiisp5Icbo/fXXaXvRyFNffUZvm4TRoYO/LUsSqqeEsQLUxvJbQKmg+kGt\nFnscgM1rS3bJgS5IF2GEOsg8ukp5ryjCWG899W+K5bZn580DLhc0n0dtoiUMn0oqrxQR+55FGb2/\n/BJYe2113NBQafTW14Gs0duGjz6yq6pc9clMxwXTZhbzjMte4ivz0LMaLpWU7f2khKHxyiuV4UJz\nEqSEETN/oZZJshry3bbfXs3+lumbHnW+OJiLkzDMd7MNshob0zJ3ffO62zCYuZmZz4daXnxzZj6D\nmU+HIpCBtSVbPFwL6blgityuzjNk9K7WuJWXMFzvpvP98cfhfLz8MnCGkBvzeEl16qT+bWoPXzw+\nKaJowgh1kP/8Z+qvTuQnjFCH75prYFNJ5RkAVLMOkyv+r7/OZxy1xemSMGzxyLZkql5l+wp5xU2f\nrryBgKyE0dpanwX1NGwOAwsXKpuLVDEffLBdxVtPt1rzG8t60tiY5k+qVzWkSqpWco3hwVWQdaNd\nlFxbohCzkJ6ESRgxEsa554ZdF2MljBdfdC8kWI0N44MPsucxI9w8Kinb5DYz3qIlDN0IilRJacRI\nGCbku7s6dp1n24qjMY21GpWU7FRkuWy5ZaVO3syDz6sKsKvWdJ7M8NLhQzoRMGeXtLERhqtsJGF8\n8w1wxx3uvBYJSRjLL58ljHvuAW6/3f1MLfmLlTBkvyQljDXXRAVknG1BGLcBeIWImonofAAvI9ln\ne0lCXj2maZx+6aXsfS26yUZy4YWVHUm1hAEoLxPfUht51Bqm6B6jv5bn8+a5l35mdhOGrIw6zfnz\nw15ooXuASnPKlGyZxxCGXtPIjNuU2GxutVLCMPP1xRfptTwSRqytxhav+czs2ZUrrvrUEY88kq6u\nrMO8/baadxEzk9o3D8NUT9oI45e/VF5pkiRsAwDX4EW2j2oW2vv00+xK0z7Y5lGYhGGbm2NrEyGV\n1Pjx7vZtkzDkfSBb9k1N/uVW2pQwmPkiAEcAmAlgBoDDmfni2pJtf5gShrk09JNPZsPpf7PSLlqU\n2hEA9weZP7/yo773nt3FMNSRfvJJeNXaGP21zOt//1vpFivz0b9/5TOAKj9zXSdt3Jw6NSXialVS\nc+Zk09Tfy9cgBwywXzeJwaYfl0Zv8137909nbLvS96mkeve2r6klEbJhnHuukhwkfJ1Ba6taa0zj\nyy+BTTZRnejKK4dnMfvmYbgI45tv0n1ZALUCq5xRrQlj1qx0RBwjYejjGFWtxu9/719VwRWPlDBs\nXmQ2wpD1OZSvgQPTFYElqiEMqZKy4eWXgYsuUsfPPAMcX8OkiJi1pFYDMA3AgwAeAjA9ubZEQRby\nxIlhJjUJIxTOlEgkxo6150OiVy/7PIWYNM3j730vdVF0Ia9Kyhe+tVV1FgMHVoazGdEmTFD/xx4L\n/OQn9rTlNR9hmMuc1KKSchFGrIQBZFfNtUGXz6JF6ZwHHTZmGW6XhDFnjorTtvfIppumYW3lYnNL\n1v+2ZTpco+eQhKG/zWmnVS6n8te/pseaMKSEEyNh2CTCUD3Iq3nQkBJGx46V3pL6f/ToSgnPJmHY\n8mkbPLhUUprYXSopH2H84Q/pop433QT85S/usCHEqKSGA3gcwGMAngYwDkDk4s9tB1kx+vdXBePC\nk0+qzXuAMGGYowe9K19MPiTyeCeEOlLmtBPIs3yFT8JwjU4B1Ykxq4pphpMNw1SlyU6jKMLweavE\nzsPQYX0d0OLF1c1g1/VpwoR0M6c8nZaLMFZYQXXEtjxp+5UrHSLghRey13Q8NkOt/j/zzGz7MBeQ\nNDvF1tY4t2ub40IMYWhXdUmMMXM6YuFSSXXs6J7Au+mmacfvkzBs+bD1US6jty8eacOwhZNlZLr7\n5kWMSmpDZt4o+a0FYCsA/6ktWQUi2pWIxhLRh0RUsQMEEQ0mollE9Gbyc66SaxakbeXNs89Wqpfm\n5tQlMlbCCEksNnG5WsR0pKE9B2zPmlu2xkoY+ripqbK8bOWnr0lDp5mf555L1WrmvZaW7DpcsRJG\nSA0Xo5LSjWv06HjCsI3IQ2Xris/smGQ6H36YxrvddpUb6vjqi54LYYb1dbhXXeXPu00lFbMXhi3t\nkEpq7tx0kyOf1ONKKwSfSspGGOZS5mZaRRm9XYM80wHDJ2G0KWGYSLZJdWxMGg8iagRwLYBdAawP\n4GAiWs8SdBQzb5b8LnTnK3tuVtznn1dSxa23ZgswtDZ/bIXT7qr1IgzTiBwyjtr0tieemB7fc09W\nRRJDGDYJw+YNFEMYP/oRMHSo/d5ll2VX+s1rw4iBy61Wl/OIEfk9wmT+Qvsyu2ZlE2WN2qZ0pb0B\nX3yxctntPAbNkFSn4/PZXGwqqZhlP+S7ymdtS+BrwnBJw7ocXVunVith6ONFi1SHrN/VlnezTsSq\npACleZAuui6VlJlWHpWUTLtWwgh6fRszvhugJu3ldGK1YisAH+llRojoHgB7A3jfzEI1kZsdqt6V\nq6EhH2HEerjo0XuRhOHrAOTKoL40XfcOPjjrgheaN+AiDNtz+lrs1plmHuXy5rXYMMaPz5KkKWH4\nvKSA2lRSvlGi7Zos8z/9KXv9qafUsbl6gC2OvLYdn2qPWa0/5UKtEoZMe8oU5U1lQs7Ml89rNZy+\nfuKJqcutLa0YuGxZss+IJQzXd+jQAfiP0M+89ZZaRuWKK9JndTxjxrhXcFh1VdVOtBrQpykpkjBi\nxml6xndXAB2hbBl715YsADVzfII4/zy5JsEAtiWi0UQ0nIgsVSoJGJAwNDE0NMTpTs14QxXPlABq\nQUyaLiOcKy4bQjaMvBKGec22ymwoH7awNgkjpmO81XD+NiUJm4QhO8F6ShgmpKrB7Iykl5MtLxp/\n+EP8bnaxhOGDzYZhtrvYfUukK7SETW24aJFyDZaQXlmhtGwgUpKthiQISRgxE0clYdjUre++m03X\njEMPtDbeuHLRUDko22QTdRySMOQ30Cs2VIuYeaXvMfN98gIR7Q/g77UljZhP+QaAAcw8l4h2g/LS\nWtsW8OabmwHoHbQGo0OHwZn7UpURWhMnk8kllDBi9x7wxZHHhrF4cTxh2OJ66CG1qdSIEeFZ6SZR\n2SSMGJWUOevVZcNw6Z6rkdxMHfOiRemsZVsctvkq1eq/58+vtFO58mqzI+SVTmJUUraF8Wz1W3vT\nmbDNMdHzPIj8Hbj5XAhSIpZqOaJKlZSPCGWaNkcASapmvj/6KGubMtd+kvVLx5PHhjF//kgAIzF0\naHX1LEbC+K3l2tn5k6rARKhlRzQGQEkZ34KZv2bmucnxEwA6ENFKtsgOP7wZQDOam5sBDM74sz/7\nbKofzUsYsRNeZAWrFXkkjBB8+clDGFqXK59xrZJpK9OnngKeflrFac5G9UkUJmHErOkj8+xKp61U\nUtOm2ddTkmFleB95xuTB5nZrg6+jje1IYlRSvpF4HmcA+a677KLSbmqyx/Hss+kk1Grbo7kulG5v\nPslM5tX3jrLt/vSn2XumpOTyfpJElE/CGAygGeeco/vKfPCtVrsbgN0B9COiq5HaElZAMTvuvQZg\nLSIaBGASgANhbP1KRKsCmMrMTERbASBmtmrG5faSQHaks/326XFeG0YREkbM8tiuNPUez2Y+YlVS\ntRCGvLZoUaWEYRJITEdgUx34CMOcQBejEtAwR7sulZRrZVnfe9hG8jKu0Axmfc/siGwqqTzISxg2\n6aqtCCPm3WxhNQF37mxXUTY3qx0vr7yyeolfqqCIKlc68Nkw5PHixZVlIs+lYwhQORA0+yc5gNVh\nQzYMW39XLZH6xqmTALwOZa94HSlhzIbaRKkmMHMLEZ0E4CkAjQCGMfP7RDQkuX8DgP0AHE9ELQDm\nAjjIGaEBl/GtXiop0ygmod0BfVi4EHj/faWXlGlebMyp16Pi2DWHirJhaMKQFdisqLKRuGBbXpk5\nSySmSG+eA3Gdms9jR0pHptRkpmXmdexYYPfd7fHquGxrSplwORrUohYzN9JyQeevFhVqrV5SMZ2W\nz+nE1VHGTKQ14XKrdUkYNsRKGD7HADMfJmHUKmHE1EsfnITBzKMBjCaiO5m5Lnt4J2qmJ4xrN4jj\n6wBE7HdWCddHufxyYHOxOHveiXsu5LFhdOxYacwbNkxN9Lr4YuDuu7Npy+O8hBGrVvHpuQFV0Rob\ns/l2NVhfmraZxcxA1672dF0SRswGSqEyshFGSMJg9q+2quPS6jBfQ3appPIQhnk/VsKwjTTzShjV\nGr1j1byAf8DmUsvKMihKwvCppMw82ghDhveRqkuiMM8lYYRsGDLtWiUMpw2DiLRR+w0iGmP83nY9\n196Q3lA2TJ/e/l5SNjLTH1BvCmPGZX7g2H0TbBXj8cezzzDby8EmYZgdjG2U7CtTF2G48hwrYbzz\njlpjS8IM47KNmBJGyB7lM7ibEkYsYbhGpKH6ZIa3EYaNADSh+WYIu6Djawsbho9cZPpFGL0lTC8p\nUyUlJVQfYTz4YGXcPgnDtLu5vo9UZy4REgYAvcHintVF3T7QI2D9MeWMYY22tmGYsFUYvcCZrDAy\nruHD7XGE8me7fswx6l9W+JBK6tln1S5k5mg8rw2jGsKIsWFstFFlvCHpSt8364Ce++GSunwjcC19\n6DgPPNAd1vZeeSWMGMKwxWGbEEcEbLNNpe+/C42Nyqi/8srAIYeoztGcUV4PG4ZMP/b5EFwqKT0P\nw2xvRMC++6pj21pS+viwwyr38PDl29Q8mB17e9swfBsoTUr+P7X9qkuu/jBdJW265nrM9PbZMEzY\nCENPqJEVRsZ19dX29ELifWsrcNJJ2Wu2kbVNPSWvvf66XcLwzcOwwbamlllmPhuGNHKaC9yZ8I3W\nZWcr8yuXZq9GwtB50l555nwBCVv55yEM2/3Yzt7Vcbz8ctzzgJLWV1kFePhhtTkVEOe9V5QNQ6Zl\nIya5HH1emG61us3qiYxEylUcqFwWyDXY0PDdi5UwQiopmUaREoZPJTWHiL52/CzjxCUDJmHY1vxv\nS7daW+XwEYb88L6KZUo0PsKQRkDzno4j5CUFhFVSrudC93xqEZ8No1cvtUS8C6F49bksc2k0dr1H\nzMhWzjCPyZ+rnhWhkrLBpZKKgV6GQ8/InzrVraZySWkLFuSzYYRUUjb06VM7YZhGb/M+YHcgqJYw\nXBKFeV6r0fvjj1OSzwOf0bur696SDL1ZjP4oNq+RarykQpXbJTrHekmEJAwTpkTjyp/PmC1F6ZAN\nAwhLGDE2DFs+fYY+c7RmqqR8UkY1EoasL3KJEolqJ9b58ueSMEIo0ugdi5NPtscFVHastvY3caJS\nwcq9OFwIeUlpVGPD6NLFPUPc5VZrS88kDJuE4bNHyms+SUHelyop10KatrzqyYmvvqpcj/Miagk3\nItqciE4lopOJaPPwE+2Ht95S//oDFEUYoZHK/fcD++1XGc4Wv40w9EeVbqe+xuxavM6Er3LmsWEA\n8RJGkYQRsmH4VDA+wpDxyvwuXJjGbdu3WY84i0CtKinXd4uBzYZRLWQcMSopLaHE2jBaWuwdu5nW\nK6+oNanM513o0iU9NgnHZfQ27wNp/R0xIn3GRxi2POm6YHb85vfR/ZkpYcTaMCQx1WWmNxGdB7Ul\n60oAegEreGioAAAgAElEQVS4mYjOzZ9U20CPOPVHsXUottGBC/KDrbqqP+wDD1RPGDoduepmjErK\nls/Q9WnTsvdiCSO0vHmMNGZ7J9+oKuQlFePiaotXNmqdflOTIv2iJIgQXCopmb45sUuiFsLQo0tz\nABCDHj2yaxLppTqAOHWdRqwN48wz1aZhJsy0tt5abdoVm0bMukr6e/hUUvo7Xn55ei2vhKEJPEQY\n8no1brUaLS3VDXxiHjkUwPeZeSgznwfgBwB+Hnim3WCKh7YOu1oJI6ZBVauSso12QyqpmFnJvjgk\nYVSrkrKpVfJKGKahzydh6I4uhjBc+l+dV7PMd9ghboRcxKgcUJtxmZssmWW+xRb+OKrV0etl0qux\nYcybly58B2RJLnbJGv1cCMxuG5xNJZXHBtTVo3SX/UjIhmEbmJjp6vIG7KpOSRjS0yyGMBob/e3A\nRgx1kzCg1nySi+IuB2PNpyUJtsXKTDzwQGV4F2RDdjFyr16V4X3x+wjDlrYNtaikzGeKsmHkSVPC\nHB1NnZqNSz6jpSOd7v77u9My82GSsYuka1nQMQ+uvjrdxlWSd55FAPMYjm2oRsJYsMAtYRTtJeVr\nd7a0HnvMvTeGiXOFnsR895EjlUt+jA0jRBijRgFDhqTnppstkLaBhQuBbt1SN3sfYUi3Wt87uwij\nXhLGbADvEtEtRHQLgHcAzCKia5I1ppYo5F3yIM9Mb1eDkittthVhxKqkYiaOFU0Ytdow5PLPpoSh\ny65PH3carrRMYjBVUrESZKhudesWjgewz7kxVVIx+akF1UpLRaikYvLuK49ajd7du/vT1m65NhvG\nO++kx7Y2JtN1OU9I6LqgF/kMuXabE/dsm09puNpcNRJGjAD5YPLTGCmOCxprFQfbBjY+5Jm45yrg\nTp3SYzPd55+vDN+WhOEzCst3s5VDXrfaomwYZh7kM/p9Yjomn73FZvTWE/Z8Dck2qjfD/+536Wqp\nPkh9tc8bKE9e8iIkYfTsqZZn15KQhjQY11slFUMYvuddiJEkXRKGhE0NJtONWd9LqqRiCEOqpFpb\nK439ErYBXLUSRvDzMvMt+aNtP+SVMEKEoT9ELGHIpT0AYK+9KsPb1pLxeTPVghg3S1PCcOnTTaN3\nvSQMM2/yGU0YMRsFmWnJZ3yEEUIoTO/e4TiAtJOYPTtbb/OM/GINxy7E2DBsJGBKGL6wGkOGKJvM\ntdeG86XhkzBCE/eA6gjjRz9SKxvo5202DImQSiqmDUqVlCQMF6RKKvQNXYRRLy+pPYnoTSL6ammY\nuFe0Surmm9P4XIwsCeOss9JjOQqTKELCMO+5Pn7M6CaPSkp2upowYueE2OIEKo3eEi4JQ+djs838\nz0qYEoZppPctCaIRM6qX9cEHLWF07w6cc04afx4bhm1iqi2cCyEJgzlMGLJT9436u3QBVhK72YTK\nUX8PV7uLUUk9/LA/fg05w12+r8tLSiK0HlcMYZgSRghz56bvH6ovbSphALgSwM8AvMPMBfmH1A8x\nRm9b+BB8k1xcHYQrD7YKYQvrewcz3w0N9ncpmjCk6KsJw1RVFSlhmDYM/T76mfffByZMqHzOlg+Z\njm158xjCAML3Y5b4BrJEqWes57FhMAPrrRcXTkJKUqFOW45kJeTe0LKcQ2SZZ9kTvUFSLTaM2LyY\n6Wq4bBgSIQkjj0qKWb1XjLpMv3+o4y/ShhHDMZ8DeHdpIAsgnT3a2gpcc004fLVrqkhojwYTrpFF\nrEoqdmQIuCuNbwaoTCeGMMw0tA3D3L6ynjaMefNUp6Qb2Pz5wI47up91pbPXXtlJTEBxhBErYUjC\n0J1KUR5YEmY5mCNoH2THJCHrsDR6+96dyL0PiQ1NTX6JK4+9xJWfULytrbUTRh6VlE4/pjNvDwkj\n5pGzADxBRL8lojOS3+n5k2obPPOM+mcGTjklHL4IwojtIDRiVVK+vJn3XJUm5v1saparrwYOPjh7\nzUYYUiUVk2Y1EoZsfC0tiqBlZ2vzQX/qqcqRnU5n3XVVnk2V1LrrpmXha0whX/9qJAytaguppPbY\nw52uCzbnBY0Yt3JbxyzrcCxhAFlHCdd3P/NM9Z9HwqgGIcLQ0q3P6L3NNuHBVh6VlE4/hjB0HQ2F\nvffeymv1lDAuADAHav5F1+S3Qv6kKkFEuxLRWCL6kIjOcoS5Ork/mog8GussijJ6x6A9CMOUHFwd\nXLWEcdNNle6ALsKw7YvgQszEPTO8+YyUMHQ+TOy6q9pDXEJ/6w4dKudhnHWWavz/+EdWajIRI4HE\n1AfTwK7JbezYrFuxiauuCsdtwqZa1KhWwqiGMIjiCGO//dJ81mrD8MH1jL6+eHHY6G1b/aAWlZSO\nM+Z9zAVIzTz4UE8bRh9m3jl/1H4QUSOAawH8GGpy4KtE9Agzvy/C7A5gTWZei4i2BnA91EzzIIq2\nYfhQK2G4XDV9ZGa6y7o+fgwhxrpnxqqkfGVaq0oKqJQwYjsLF2Ho3QQlGhvdRBZyOIiRMFyEAQDP\nPed/TiNm8T4gXiXlMnrbCMNmGAbiJAxdPxYvBjbYQC2IJ/etkS66sV5S1SBEGC0tYaO3bR0nkzBs\ne8CYMLc+zkMY1aClJV4SlojhmOFEtEv+qIPYCsBHyf4aiwDcA7V/uMReUOtYgZlfBtCDiAIrOinY\nOiZbxa+nDcMFkzBso2iXEVvDVMMUrZKyxWdWMFPCiCGMvBKGLW8xEoYNOl+aMGR+zaVWfBJGSCUV\nM4Awn5k3L58aIg9s7tEaMRJGjEpKI69Kqm9f9yxqLWEUbfQ+7ri4Z6SE4VJJ2dqpSRi2TdxMyPrc\npUv8JFL5rzF+fPjZetowToCyYcwv2K22HwDp3/J5ci0Upn9M5LaGYGPUJUElVQRh1KqSinHTNYmx\nGpVUERJGp06V7r0xCEkYMh5dno88UhmPmR9m5Rqr11eqZuSrt8ANoZpRpU8lZc6rsT0bUklJKcA3\neLKppGz1Vl/TRu8YlVQIhx+eHh9xRJofH2IljNDEvRdfDOdPqpi7d89HGCYGDqy8ts02wPbbp+d1\nm+ldx30xYv1BzNdyPNcsjgeDeXBlRJYCcnWoxx4LfPCBWlPGBe3xkVe0MwlDj2TMuH2zN4smDHPG\ndqxIXG8Jo1obhg0+G4YpYejy3Gkne54kdAeh37/a5c9jnmsPCSOkkspr9Nb1Q3da5veLVUnJdhSq\nA7femh6HvIukSiokYcSopGIGpZIwevSoTcKw4cADlQpz1KiRAEbi7bfd88R8iBoLEdGKANaCMnwn\nmeRn8yeXwUQAA8T5AFQuamiG6Z9cs6A5c2YrRNtHcH3Ms88GLrvMTxh6dJFXwjAJ5r33gKOPzl5r\naPB7ecUSRrU2DNdo0wxj2jCY45Yjic2jrUPr0kXtG9Khg59sXOnokatPJaXL07aDnPkOmtRqJYx6\nSRg+L6lqbRiyA73//nQin67bN96oZnSbW9TaJAwz3VgJQ7ajPO7IeQgj5CUVQxgxqIYw8qjTU3fm\nwQAGY+21gf79gddeOz9XPmNmeh8D4FkA/wRwPoCnYPbO1eE1AGsR0SAi6gjgQACmAuARAIcl+fgB\ngJnM7Bl3p4glDFehy4rtgq7IeW0Ypnh7//3uuF2ItWEUQRhyuWW50KJNJTVsmH92bTVutSYGDQLe\nfDMli1okDK36M1VSvpnLNpWUvLbKKnH5MRFDGNWQUSxh2BBDGFOmZO1DgL8zDhGGLPsZM9x51Gn9\n6lf+dzCRZ/6CJqwddnDHVQRhyPaclzBsadkm9spyrKcN41QoA/WnzLwDgM0ARC4g7AYztwA4CYqA\n3gNwLzO/T0RDiGhIEmY4gHFE9BGAG6DsKVGwjfzzEkaIwWPFcBMunf9xx6Vr/4c6D1nBOnVyf/w5\nc+Ly5CMMuXGU9pHXYUzC0I1YLvkuUYSEYW6mUwth6AZvdlwuwrCRq5Yu9Lt16QJsvHFcntZeOz2O\nacBF2DBMdZIvbpdKyBxxa6lSh5XzXCRMlZSPMDp0AC66CLjnnsp4gFTCWG65fKPtWMKQEsbmmwNP\nPmmPy0zbZncLQUoYK0ROWvC9syl5m4SxaFH95mHMZ+Z5AEBEyzHzWADr5E+qEsz8BDOvw8xrMvMl\nybUbmPkGEeak5P4mzPxGfNyqMcr9EmIIQ4eJkTD0/byEYXYMOp6GBjW6sIWROOigrNpn+HB3+BiX\nvpCEIdeK0o1tn33cKqnDDnMv8e2yYbgqry28ufNhXi+pjh3TvGoJw9ZxueK2SRhmJxHKk77/wQfp\ntbYyemvC6NYtTsKw1S1JOo2N6R4msv3YVJimhGGzYej0fLOr5X2b4dmHPIQh399loC+CMOQAsEOH\n2iWMEGHUU8KYkNgwHgIwgogeAfBp/qTaFrpDkDuWxdgw5OzJ0KhFf5S8KikzHzqdhoY0fV/a66yT\nrRANDdlNhyTkJjEuhAhDdmT6+Pbb08Zv7mDo8yN3qaRcnYOtMch1jMy8+qDLtHv3rIRhm4fha/Ah\nlVRMnmyNta1VUocfXj1hyO8lV+eV7SckYbiIWkoYPmgJwzZ5zodqbBiu8PWQMGIm7p1xhl8yNyf2\nmgPguhEGM/+Mmb9i5mYA5wL4G4B98ifVttCNOOSrbdP16bAxyyYAtUsYkjB0Hr/+2v28NgZq1DKB\nB7AThhSLbYSh0zS9pEKE4VJJudwWzeVJgDhVowtnnKFmgYdUUr4GH1JJxeTJRg5tJWFIg36MW62r\no9SQ0qRUSZmIdauNJQx9P7SntYlqJYxqCCO2U85LGKuskpZjjIRhW2KnXiqpb8HMI5n5EWaOWNKu\nfaE9V2TFdm1VKJFHJaXhIowLLrBf9xFGTAUzO1fXM+ZI3AUbYciZxzaVlGz8eQgjr4Rhg20CYSzW\nXz/Vr+vOwKZLd317mzShy+/Xv05tPKE82QiyrSQMqWKx7YMi4ZIw5PvJpc59Kil93SQMM37pJeWD\ntGHELLKpkcfoLSVol1RoEsYdd6TtJ3ZOjlRJxcz0Ds3TMglj4UJlh9Gop0pqqYOurLVKGD7CkJ2W\nizBkQzLzJyFtGDGdX+yHjjWe2Qjjc+HgLNOTSyqbRu/WVnXuGyHZOpGxY/PNZalFotKN8eGHlQeO\nbvBmI3VJGD6j9//7f8Af/xiXx7aUMFzqxqIJ47HHwiopmWZoHkashLHccvYFKF0ILdqn833EEXEq\nKbN833/f745vw4IFWaKMqT++vX9shCH36mkTCWNpgm7YeWfdytGHj8Gl3cLV2bk6gIYG4OKL0/O8\nhGG+k4vYtAE9BLMTXLwYmDQpm18NU8L41a+At95S1xYtss+aDuX1scf8KrhQHD4CvdrYdV7m7dln\n86ukbPppKa1oVCNh1GvinvmMJH1Zx11qnZAzgCaMzp3jjN6mDcMVd6wNY7nl0lWqYxCSMGS+Qyop\n10jf5SXowsKFWSN+DGHkMXqbhDp5cilhZCB11BoxBRRLGFLd4yIM3Slst132ekMDsNpq6XkeldQh\nh6hp/hIuwlhxRX9cGiZhfPKJmtQj86vDycZjLts8fz5w2mn5bRg6vhD0ZkEuNaINP/5x9tw2Qc9G\ncnkJw3Q/XZIkDLPjlTYMWZauSZCxEoYsW1c9tqmkXF5SsSqpat3aY+qo/K6ucrD1E3k744ULs0b8\nWgnDVNHZCKOadfScr0VEc8TaUeZvid2iVaOaUR8Qp5Lab78sYbgaur5uM9KaI3p9PZTHO+6onBjm\nymfs1H+bmmXnne26Wzk6/fLL7DNab5tXwgDiJEGtg9XlFfOMbRVamwdYrIThUkkB2ef33NOfr7a0\nYbgW99NkqeEijFgJo7HRr5IyJYxFi2rzkpIqqTwIEYZsN9VKGHm3il6wIGvEz0MYNthUUkXAWf2Y\nuSszr+D4ObzslwzoDtmUMGIIQ1Ymn+FTkoCrEbsIw/SJnjUrve7rEP7xD3t6tRLG3ntXxrHKKmn+\nXV5SNgkDcHvW6Hs2xHT+UpUBpCMy33e1qWNqIQybhGEjjHPP9Xd4NnJoq4l7LqO3rVPxSQoakjB8\nKildL3wShowjViVVtIRx3HHKbXyzzcI2DBthbLFF/tG7lDDyGr2rUUm5ngsharxCRJsQ0clEdBIR\nbZI/mbZHrNHbhO64fCop5vTj6kW8bKMznbapsjIJ44EH0uu+PP7sZ2k4CRnXZ5+lx7aF82x4663K\nytOjR6WEYXpJmYghDBe5xXSWJmHojiIPYdhUUkD8PAzbYogxWwGbsBHkf/8bfq4awnDNNTI7O1sd\ndkm98pqWtmXZ+ryk9HWb0bupKc1fvSQMHwEAqr1uuGE6OPB5Sdn6iZAHkw3ShhErOdcqYeSZ7KgR\nbKZEdCqAOwH0ArAqgDuIKGLz0/aDrpSm0TuvhNHVs06vJoGNNlL/volYNpWUy4MkTx415Ifv2zc9\nzrMapVl5evTIjr7NtG351KOYaggj5r2rIYyQSkp2BtVKGC7kyVcsqlFJuTwBTZWU7LQ0YgjDppLS\nrsu+uGzOBlIdk8fonQchwiBKB3XVSBimbSgGCxbkkzCqNXrPnq1WinA9F0JM9TsawNbMfB4znwu1\n490x+ZNqO0jCyCthSH3ldddV3u/XT60r7/Idl/CppFyGsmo8ZWQnbNpsYkekMYQR0udqCcO1/pCO\nw4ZqCCPGFbdalZSL2GxL0cemLaEHMq++GheXRpEShtmxMbvtHb582IzePglD5sssd7ksRogwdBnW\nizDk4oOu8CEJI7auVGP09pGSS/pYYYWsG3xexI5XWh3HSyS0ntRUSeXpjInUEhJycThAzU845ZR8\n6q2QSkpez0NqGkXM+jbz0727XRT3lWEtKqmYfGvPLd0YhgxJ03Mhj0oqRsJYuDC+ofmkCH0vr9t3\nNRKGTyUl31OqgzRiJAzbyNhn9Jb5MtObOTPeS0qXYT0IQ8+vkIMflyt0SCUVQxpSwijCS0pixRXT\n+UES9ZIwbgbwMhE1E9H5AP4D4Kb8SbUtqpUwTJVLzAgrFFcelVQe11+N1la1rLgtX9VKGHIF3FhP\nMy32tra695t2VdKY9z7sMDVfQ3eAZ5+dpudC0V5SCxYUK2HkVU1VMyDwqaQkttwynjBkOP0usv7G\nGMu1DUOHXWcdtZ95rIRRT8LQgzopYdjyY3OOqVYllWceRh47yZZb2l3sC5cwiKgBwMsAjgDwFYDp\nAA5n5ivyJ9W2KMKtNiacD7FeUvJ6NRKGTwUUCzM/TU2VXlI+yUHG4etQL7nEfj32vbt2te8/4HtG\nwkUY5nVXY1q4sBjCaEsJwzfTW8Zrm0gXM2CySWmxKin5zFprqb1XfIQh5we5nEpCyKOSkm3Llk6M\nhBGDoudh+KDDFy5hMHMrgOuY+XVmvoqZr2bmN/Mn0/YowugN+PfLMPGXv9jjMgnDtpwAUBth1Aoz\njg4d0vyfeqo77c6dgT//OXutmooY896ub5JXJWVTscVKGBMnpt5qIfg691//Gvjtb9Myjt10SeZx\n0KDK+3ffnT0/7TR3vmyEYZNQYzovHUcelZSeh2GLB7CT6YYb2sObEzR9yKOSCkkYpuOAvlaNDaPa\nmd4aLieXvPZEH2LGK08T0X5EtY5h2xZFqaRchGHrDGyjWaCSMPRexSZijdRmOpttFh4J7rGHP06b\nhKHT6devMj6Z9s9/7o8rBrUQRluqpKTbcgg+wth5Z7U8jE4/dvFFmcfnnkvJvKEBOOqo1ANGp3Hi\nie58mfY9l4QRSxihiXtm/hctyn4Ps+25VEBmukSqHsbCJmGZkF5SOpxNwihSJSVd+mPK/JJLlGOO\nLuc8dROoH2EcB+A+AAuLmulNRCsR0Qgi+oCI/klE1lWPiOhTInqbiN4kolfi47cbvestYcgKeNFF\naQUwK75LwvCpfGwNRWP11cPPPPZYeryPZXF6n0pK5i/GO6wq3WhETXSNDH3ifx6VVAw++aTy2hNP\n2MP64jRH0dUQhtT/b7YZ8Le/qeMRI5T31eOP+ycImhKGVr/IhfNiRqc2lZTte9okDNszPsKwSSSy\nHGIQkjC0FLp4sWrHep0zV36KUknlUZ83NAA/+hFwwgnZazbkdUDxphsKkMz4bmDmDgXO9P4NgBHM\nvDaAZ5Jza/IABjPzZsy8VWzkenTT1hKGTOvss916apeE4WP8kLeXy5MlRGwaPpWU71nbSLJWlZSL\nADVOOAEYNSo9r4UwXCqpPHDp0EMjU50fXxw+EKVqKZn3H/9YGTo7dLDb5FzODFrC2H779HrIbiXf\nwVRJbbONskuYedZwSRghFZCZrlYhxSJEGJ07Zwd1776r/l3fMdatds89gWOPtaeZlzDk++o0Ygab\nEnWRMIioYh1I27Wc2AvArcnxrfBvyFRVM652prdJGHlW8HSppMyK31aE4YItrLn3t83FEgB69gT2\n2sufl1oJI9T4l19eja40fOK/GZdJDC6VVGx5HnKIe2kKHbdtXxJTwqiGMBoagJNOAtZcE9htN3cY\nmY685pIwJPIQhqmSuvNOtdS3LT9AasO47bbsPZ1eaEXfIgjD/HbduimSk44puvOPlTBcKqnOnd1b\nDsjlzWUeXciztIztmwIFSxhE1JmIegLolaiQ9G8QgH6u5yKxKjNPSY6nQM0gt4GhbCivEVH0ZMH2\nkjBchGFWfJ9KygVbQ5GolTAeeih7biMMPbHr4Yft+ZLh8sI2Ajbf0/WORUgY5vXY8rzjDjdh6Lht\n+6IUJWE0NAAffgj87nf+PMjOLiRhSPgWYbTFJzvjhgb/N9QSxvrrq3Oz7YXaWBGEMX06cOGF6b0B\nA9J0YiaKhozeEr52YUoYIdQiYdTLS2oIgNcArAPgdfF7BMC1oYgTG8UYy28vGY6ZGYoYbNiOmTcD\nsBuAE4nof3xpHnqo+u/TJ37HPRNF2TBkXD4JQ04M9G14ZPPsCeXHdX2lldzpaLhUUiYaG5Wrq9af\nA7XbMHyqHBvqoZLKQ8Ch/VBiJIw8Ow6accSEsbWDekgYIYOyjGvBgmy+YuZAuYze1dowunRxSzJ6\nXpGekGr7RjZpohrCMMsiBJ9DQGzdraadOr3AmflKAFcS0SnMfLUrnOf5nV33iGgKEfVm5i+IqA+A\nqY44Jif/04joQQBbAXjOHmszxo9XG/pMnToYzIOrUkmZYnE1hHHDDepfV0Tblqr6Y8nRZ8+e7nzl\nJT5XHk87TRm9b7zR/6zN6G2D7miPOgo4+mh1rVYvKVenU4SE4VJJxbgyutDPIW/rtG2EUYSEEQOb\nSspMW19zSRjdAhZLSUA2m4lMQ16fPdtOGLr+2BZDrIcNw9aRNzamalr9b2t3TU2V+ZQkIuP2DQoW\nLgy/wwEHAPfdl6arodPQ+evUKSU5IFvmI0eOxJgxIwEAb7zhT8+GGKP31US0LRH9PyI6TP/yJ5XB\nIwB+kRz/AsBDZgAiWp6IVkiOuwD4CYAx7iibscMOzbj00masscbgb1VSeTuBImwYgwdn47JtYGPr\nVH0j/03EGsGhkZfv+qBBcYTj8pKKSbcolVQs+vRx3zPfoaUlm79abRiAInpXpwP4JQyTMF57LT7d\nGNhUUjZC9kkYu+wCfPRR5XUNG+nG1MdZs+ySpW4bcl952/MxhKHtI7Y4fHVXtlHTvidhIwyXhNG9\nu7ttzJ2rHBWGD3enJeuRTSrS7+Nb7n3w4MG4665mAM3YdNNmd0AHYozedwD4E4AfAvi++NWC3wPY\nmYg+ALBjcg4i6ktEjydhegN4jojegppt/hgz/zMmcm3DqGamdy0qKdPLI8aGISuQnMVq4p+BN4/t\n4GLDSQ8bH+phw8ijXrvvPmBVlwXM8syMGW79eyitvNBx7rFH5WDAHOXrOuJbIdlETDnb6qCLMFwS\nBhHwve+F05Ck61KbmoRhkzD0szbCkJDSjIswpHOEzIf8dy3To6GX3bGhqaly6XAXYUhJ7bPPgN//\nPnu/Y0e38wKgBnv335+mq6Hz73KgMOvyxhtnn8uDmLHcFlC2hBOY+WT9y59UCmaewcw/Zua1mfkn\nzDwzuT6JmfdIjscx86bJb0NmdiwqkUJWVi1hhHT/JmJVUj6DnEkYNgnj2GOV6kp/tMmT06XSbQit\nl1M0YTQ1qRnccoQWalgaJ50Ul4YrnjwShqtx2nD55crVU0p3ck2nmEUct9giX94AoLm5cq8LU1Uk\n12MqEjaVlHkPUO/rkjBscM3D8KmkzOsmYehne/YExo6t3KALyHaGMTaMmCX+bd9d5ut/PJZTmw2j\nsTFMGKutphYxlQgtE8Oc7jrp27ExdkOpeq1W+w4Aj9C/5ECOHMxZmvK+Dy4Jw9xy06eSco0eZRq9\neyvS0JW1d+9w3mzw5cNUs+y/P7DrrnHxNjWpvcjNWdyutDR69AC2ipwxc++96bFNwogZAZleKocf\n7g572mlqBC/j1UQc21E/+WRcuFCc+h2bmrLu3+1FGD4JI4SQSurgg9Njed1sn7KTXmcdt0pKu3Xb\npBM5ar/+ers60ETsQMiGkEpKxt29ezacKRXZPJ9MyHpjhnV9a1e/Vy8JoxeA95IZ2Y8mv0fyJ1V/\nyMpalErK5bIWI2G4PGDks9V8NBtsHe6BB2bD3Hef8tmPQeyieGY55JnhussuwMorq+MYwnDZbmSa\nN98cTleOrPRoTEulvrSANL8xkGXTrRuwwQb2e0C6tEUewoglVCC8h3hDAzB1KjBuXHz6ZjymW63G\nXXel18xy9TlzhFRSNhvGfvtlw8S0eVmOZucbgoswTKlj9dVT26aGSRgxbc5GGCZi3dHrRRjNUBPr\nLoKyZVyW/JY4mIRhitjVeEkB8fMeNDHE2DA06kEY0lsij4E8FCZmJCY7b+1f74I0LP7pT/nyJ+PI\nuwyDTcKIbWR5IOPs2BF45x13/EOHVi7GF0KMe3QeCcOGalVSMW61+hkN03nBppKScdhUUrEbQIWg\n4zON5jcZGzu4iNisk+PGZRdOlGn44pKQLs42CUOmL+Eqg3otDTISwFgA3QCsAOA9Zh6VP6n6Q9ow\n3h2WlIcAACAASURBVHmnUsL4y1/S0Y6J449X/716Vd6zuVzaGoS5ZaiLMIpWO8g0zWMfXJO9qk0X\nyDaUGFdHHf4HP6i8FyNhNDQAL78cTscVr0vC8CF2zkSMSkqeu2bX23DMMXH5yOMlZUMelZRLwpDp\nughj0iTgzDOz96rxkjI70rwShhm3WcahzhmIH8T4yLMalZQrT20qYRDRAVBeSvsDOADAK0S0f/6k\n6g9ZWbV7miysbbbJ6lNNjB8PnH9+5XUbYdg+gux89HOAffFBDddHe/ppdz5D+fjhD1PjrK/BnHtu\nXBo++CSMEGFoKdCMJ48No5oRpE0l5VMbmohV19neX69L5IIrXVOdEbsaqm+eR1Gq0WrdamUe+vSp\nLNcLLlCDvEcftc9RChHGjBlxBOxzFw8RhkvCyLtabUyeVlmlOsJwoV5G73MAfJ+ZD2Pmw6Bcagvo\nauqHPG6ZEgMG2L2RYiUM3ShjjN4aroaax71SpgkATz0F/Oc/+Z6PQYxKSjYUfc90FdS2FebUJTHG\nPddn4M8Dm0qqU6c4GwZQvX0H8LtOu54BgHXXzZ7bJrX54rNNvpNp5e3gbGUlVUOxbdDXuX3/+2ob\n3p/+1L73R4gwJk6s3oZh85Yyw5rpafgkDB8x+zrwo49WWpAYG0asi3i9bBgEYJo4n55cW+IQGt3E\nPGtrsN27V3bgMRKGafQ+6qjKNIqyYcgK1NSUfyc3QG2BCgBrrBH/jHwX0wtE2lIktO69oSFdfsHn\nnKCXfAmlHwubW22nTtnrRRDGeuvlz1vs+8R28Po95LIzus7JznDGjLj4fGlIlVMRhCGh25BtwqDL\nhhGrkvIh9D3yEoav7vj6gt69swPXIiSMehHGkwCeIqLDiegIAMMBOHYAaF+ERjd54pDo2RO45hrg\nvffSa3lUUvrjnnde9joQZ1SMwZ572tVYecpCd0Iffxz/jCyvgQOz90yPM/nM2LFpJ9ali9/Ifvvt\n7vSr+da20XGefaFjCeOGG9TyFxKhxmyrw6Z0AcRLGDoeOeCxDY5Msg8hpGN3ubPWShgyDpuEof8/\n+kg5UlRrw3Dlz9U59+2bveYidF89i1ERLfFeUsz8KwA3ANgYwEYAbmDmX+dPqv6oRcLQcBFG9+7Z\nEWMelZSu7D4pplY0NQE77VRbHKF8hZYG+cUv0lmkgL0z+M1vgEsvVX72Gl262CWMM88EHnzQnlZs\nnm2wSRIxKqkVV1T/sYTRsWPlgpJdutg3YdIw32fQILVEuFn2sYShIdcss9XDdddVExs1nnsuX/wa\nOk7XCr3m+8V+P9m2TAlJEoZ+tz593AMREzaVlEaI0HRd2FmsnOeTMGyEoecO+fKhscR6SRHRWkT0\nQ5UhfoCZT2fm0wFMIyLPQgHtD7nwVl7YjIM2vXM1EoZNL+qTMD7/vDJ8XuRxq7Vt5xmCbPCnnw6M\nHl15T6Y3YEBlZ2IShkbv3vbdAV3px8JGDKbazCyjiy9Ot8CsRt0nYdPHa5jqTNe8mbyE4VvPClC2\nJNmZ/fCH7nrnWyZEl1vs/tKxdfu44yqv2QijoUHVQV3HapUwQhMZpUpTPpOHMHQcRUkY7WXDuBKA\nbSvW2cm9JQ66YOTaQrHShqzoU6ak199+G7jWspi7z63W5SVlG9m5PlqPHukqqL5GVcR8AY2YWbEm\nfB22TSVlexeXhGG6ONvUJrHvL/f7sDXMjh39Esaqq6bSQuxs+Wqgy8fs8Mx6ktdILTuq1VZT/7LM\nFy1SzghSynB927feCo+AbRKGvK8RSxjmhDz5rHZJ1sdSyo2pH76Vgs0ysG1lDGQJI69KSvcP1RJG\nNTaMCy9UK3vnhY8wVmXmt82LybXV8ydVf+jCPPLI2lxGV1klPd5oI7uHia0ihrykYglj9uys+qsW\nCaNI5F1CwSZh2MK7VAeS+BcssC8j7lMNSsj45Xvokbq8by79/uGHSt0m74e8naqFfh9N3q4OrxYJ\nY999gW++qSSMlVZSy6douKQEnwefjtP2rO1daqnbsn5tu232Wh6cfrpaN82GUP70/cmTs/mqRsLI\nM+LP4/ZtK/f/+z+19E9e+Iq3h+deDhNh20EWjF7CYe21ld68aMRIGObCcrbRgU19Zuq9i1ZJFRFW\nI0bC8F07/XT1fWwShiTu0L7ZEqH9G2TD1F5a8voxx2TLYs016zMT3AbTBmCTMIYNA/7wh3zxyg58\nxRVV/PKdbASkbTY+mOuG6ZG163uZo+ha7HlSwujdW+2eV426cPnl04m7Zoce61Yrl1TJSxg2DzCX\nS7mO1zUAsuW5yLrq+1yvEVHFluXJVqmvF5eF4iALRjeQxsbs6LAonH02cN112WvmSMEmYTz/fFYv\nPXNmOK16ShjPP+9fHlwir4Rhy/eECdnzyy5TdgqzUj/4YHhpEcDeGFzqEA3daX31VXZp6th5GDH3\nq4WLMCSOPDJ+gUcNrc574YXU4cC0YZjwLT2iy8qcZe9azgPIGqw18tZt6bprqjxjlkqRMMO//z7w\n+OPZa2b93nvvrDegqZL6yU/yq6RsNgxdTuYcphhVpM7ziy+Gw+aFj49/CeBBIjoEKUFsAaATgJ8V\nn5Xa0KtXdnmJavTxebD++pUdmqlyMm0ZRJVi4Jw5/s5nhRX8RtIePjkQ4Y5tu+3yuZSayKOS2ntv\n9/LgMp9EYWO3hhzJ6r0Cpk/3x7/LLoqoevTIjqxlZ3bJJfZ4bPEVCVMlVRS01CXjldJVtRIGkC2L\nUIdmShi1DIZcbtsxsA1+bO7LZv769gUuuiidG2RKEuedB4wYkd+G8cMfAjvumF7T7tjDh2ffL2br\nV51nvT95kXXVt0XrF0S0LYAdAGwIgKE2MfpXcckXh6nGJq+hUaaJogp1wYJKcVx/UJ+B24Vx49xi\ndmOjMkDWiph37927cl+Ahx9Wu4S5YDa2hyr2VawdcgXYffdV/337KruDC926KVUYkN2DRH6fIUP8\n6bYVYRSVjiYM2WFL6dZGGL4R+5Zbqo2hbNd9I1uzDeRVSZlS4Mcf1+9bAHZCk2TwzTfqX6uWW1pU\nWbsM2IMHA1ddlb3W1FTpxvzll9lz6dCzww5xeTadbYqAV+PHzAzgX8lvqYLU2dazQpmw6W51JbdV\noqYmv+rFtZz2EUeoEWA1xlezPGLKRxr1NPS+BC5Uo5/+/vf9632Z6NJFbV8r3XlffRW45Rbgl79M\nr7neccgQ/8ZVLixthKHtYrIOfv21+t9ii8pJl4BfwujVC3jsscrrRGrNNheKkDAkaeRZlSAvhg+3\nD4gkYejtWzVhaOcBF5ZfPv/GSRKdOwP/MnpjFwmbTjhFoA5TyMIgov2J6F0iWkxEm3vC7UpEY4no\nQyI6K08aUsLweR88/7zSE/6sjko23UhshPH55/59fF246SalVjExcuTIzHmvXpUkJmelSmywAXD3\n3fnz4oJp7I/Bo48Cp56aLx1z3490KZeRwWc7dAC2314d5/FSqRdhmDOlXS7GeTF16khccUV2b/iD\nDlJb/776aqrOk4i1CcSWW+fOtRNGEeVuthEXdtvNPujR7efll9Md8PTKui0t8ao8DRthSDX0Ntso\nNaoL222X3YDNtehpEWgXwgAwBsoO8qwrABE1ArgWwK4A1gdwMBFFr84jPWxc2HFHVdjDh2f1h0VD\nN34bYay6atirJw/MxvDSS5Wb4gwc6F7S+aCDistLrZPqYmHrdFQ8I71hakHIdlQttMeP1nXrPSKG\nDk0ncuYFM/D22yPxy19mBw+dO6sZytKQLBHb8a2zTtgO9tZbaute/X3HjFH/eb/LWmvVThqxhOHC\nbrsp4/5WW6UOBFIlpdv7uutWLtlug40wRo1Ky+jFF1O3YRsuvjirAdCD5TZXSdULzDwWAMj/JlsB\n+IiZP03C3gNgbwDvx6Sx1lrpBDxzdLbqqmqdn733zpfvatDSkjaKaqbi14oY4txtNzXKLBoNDUp8\nDu2cVitsnY6cK/Dss/4Gp5GHrJ54IlVHFI3Jk9OO5uqr1X/HjvZ5KPWEKWFMnFi5yRGgjKuhb6wl\nm/79VdvUnnmxO0ACqpPu1EkZnR9p5z0/TaeEefPUrombbKK87wDldRUD27wWPbkyL954Q/V3d91V\nHwmjXQgjEv0ASCfMzwFsnScC3Vn27JntDL74oua8RUN3Zscf71YFtTdM9+CisM46qYFu993jnqnG\nY8tmgzjggNT4ahrrXcgzEltllTgyrhYXX6xm4tbb28+HAw7IDnKKqL/dugEffKDa4zPP+D0ATeiy\nOOccNfEsD849V+2vUY29KoQTT1Q2oksuUecDBsQPPl5/vdg8bbZZapsqYjFWE8RFra9tRkw0AkBv\ny62zmfnRJMy/AZzBzG9Ynt8XwK7MfExyfiiArZn5ZEvY+rxEiRIlSnzHwczRlFI3CYOZdw6H8mIi\ngAHifACUlGFLqw39oEqUKFFi2UR7Gb0lXJ39awDWIqJBRNQRwIEA2llzWaJEiRLLLtrLrfZnRDQB\nwA8APE5ETyTX+xLR4wDAzC0ATgLwFID3ANzLzJFmpBIlSpQoUTTqZsMoUaJEiRLfLSwJKqmqUcvE\nvu8SiGgAEf07mQz5DhGdEn7quw0iaiSiN4no0fbOS3uCiHoQ0f1E9D4RvUdEPwg/9d0EEf02aSNj\niOguIuoUfuq7ASK6iYimENEYcW0lIhpBRB8Q0T+JKDi7aKkljFon9n3HsAjAacy8AZSa78T2Kgsi\nWo2IvqbAJBvP818T0aACsnIqlCqTiegWIrqggDjbHER0OBFVuWEqAOAqAMOZeT2obZarVusSUTMR\n3Z4c1/SdLXEPIqJWIqpLn5TUqWMAbM7MGwFoBFDgNNUlHjdD9ZUSvwEwgpnXBvBMcu7FUksYEBP7\nmHkRAD2xb5kDM3/BzG8lx3OgOoUor/mkQxpDRN8Q0WQi+jMRWfa2cz7/KRF9O0+emccz8wpcpa4z\nefbTap4VeeoPYHcAf4NyquDk5wrfh4iGEdEkIpqdjMabiSjnEpYV8da1ExTpdCWiOUQ03LjeHcD/\nMPNNgLILMvOsGpL6tgzN70xEI4noqBrirjdmQw2slieiJgDLQ3liLhNg5ucAfGVc3gvArcnxrQCC\na0QvzYRhm9jXxnNhlzwkI6nNALzsDwkQ0RkAfg/gDADdoKSTgQBGEFHsPFGG29OtvXAFgF8BkHPr\nrXkkopUAvAS1bP8PmLkbgJ0BdAdQ1N71zvJJJOVasS+A8QAGE5Hc3WR1ANOI6GYieoOI/lojCfq+\n8xJtDGXmGQAugyqnSQBmMvPT7ZurdseqzKw3pJ4CILgzztJMGEt0BW0PEFFXAPcDODWRNHxhuwFo\nBnASM/+TmRcz82cADgAwCMChSbjmRAd+TzL6fp2INk7u3Q5gNQCPJuqJM81RdTLyvICIXkjCPEJE\nKxPRnUQ0i4heIaKBIl+tRLRGcrx7onOeTUSfJwSnw/2UiN4ioq+SuDfS16GIYhiA4QA2h3+HyNMB\nzGLmQ5l5PAAw8+fMfBozj0ni3JaIXiWimUl+v12PNXm/3xHR80k+nyKinsltvVbazOTeDxKJ7gUi\nupyIvgQwlIi6EdFtRDQ1kdj+L6eq5xdQ0tQL+rslaIKSxGclx78A8AIRrUpETyTlP0LrrsW3O4aI\nJiYS1xmwQIRtJKKLAPwPgGuTb3y1TbqSUkjy3J+IaBoRfQxgDyP+7kLq+zypQ7pOrUlEo5LvMY3U\nskFeENH3oPb4GQQlfXcltddPCXy7Mnm4T2XmpfIHNRp+Upz/FsBZ7Z2vdiyPDlAuyL+MDL8rlIje\nYLl3C4C7kuNmAAsB/C+U3vcMAOMANCb3PwGwo3h2EFSH3ZCcjwTwAdRotxuAdwF8CGDHJL5bAdwk\nnm8FsEZyPBnAdslxdwCbJcebQY2Ivg816j0syUcHKImpBcD05Pn5ABYD+J2jHP4DYKinnFaCEuUP\ngRpgHQRgBoAVxft9CGBNKGL6N4BLknsDZVkk1w5Pyv3EJL7lANwG4EEAXZJn/gvgSBH+OU/+Bibv\n2x9KRz9a3OudpPUigF5QKtsFAN4AsAmUVPUMgPOMb3cngM5Q++BMBbCTqAu3O77zv3WebffNMACO\ng1Kd9gOwYnJvsYjvQQDXJ/noBSUxH5vcuxvAb5PjjgC2jajvBwL4mzj/OYDr2rvdtuUv+SZjxPlY\nAL2T4z4AxobiWJoljHJiX4JkNDoMwHvMfGXkYysD+JKZbUsifpHc13iNmf/BzIsBXA7VycV62zCA\nm5n5E2aeDeAJAB8w87+S+P4ORQA2LASwARF1Y+ZZzPxmcv1YADcw86uscBtUR7gNlFQxhZl7QnXu\nI6BIwYWVoIjFhT0A/JeZ72TmVma+B6qh6d1A9Pt9xMzzAdwHYNPknktKmMTM1yVlvwiq7v6Wmb9h\nJeVdBtWhxeDnAF5h5s8B/APA+kS0KaBsW1Cd8N+ZeRpUOX8C4CVmHs3MC6A6ZrP8z2fmecz8DpSx\nNHaHkjxS0QEArmDmicz8FYCL9fOJWm03KEeOeUner0RqpF4IYBAR9WPmhcwcsxnpWAA/IKLOSXv5\nMZRTxLKMR6CkTiT/wS3OllrC4HJin8R2UKqIHUi5kr5JRKZHhIkvAaxMdoNsHwDTxPm3S7KwGo58\njkijeoIp4ng+1KhVnlvW6wSgdPO7A/g0UWdokhoI4IxEHfUVEX0FNcLuk+RLGjMZwGdwd2bTA+/S\nF0rvLfGZ8YxcznKe5300pO1tZSjJ6DNxbTzi7XGHQZEumHk6lMQjd7GfDmAIEY2G8pJ6C5Xfw8yv\nzN94xH/rPGriPpZ0NAZClclk8X3/AiVpAMCvob7nK6TcyI8IZox5NJQk9xqAt5PLN+bI71INIrob\nStJch4gmJGX2ewA7E9EHUBL/70PxLMmr1QbBzE9AjViXaTDz88hP/i9Bjcr3RdLhAN/aQXaFUvFp\nDBD3G6A650k6+bzZjQ7I/BqAfUgZhk+GGr2vBtW5XMTMF5vPENH2SDpbZh4FYBQRvQDgI0cyTwP4\nGRGdn5ChiYlQ6jiJgYird653lde/hJIyBiF1eV0NjnXTJEhtobwmgHOI6NfJ5RUAbExEZyQSzEIA\nJ3CytTIpu1NIElgNSi2mj2O8icx3TTYvxfIAtD1NLkY6OYlbpqkxAapu9rRJwKwMtccCABFtB+Bp\nIhrFzOPMsMZzlwK4NPAe30kws0tK/HGeeJZaCaNEbWDlXnk+gGuIaBci6kDKw+o+qAZ7uwi+Banl\nXJqgDIfzkap5piDsTUSOY/cDKj+HEFH3RHX1NZR6BQD+CuA4ItqKFLoQ0R4J2b0IoIWITkni+F8o\nW4cLl0PZVm4lotWStPsR0WWkDOnDAaxNRAcTURMRHQhgXQByg1LXO02D0uM7yyd5t/sAXETKPXYg\ngNMA3OErnwS/APBPAOtB2SQ2gbI7dIaSzKrFOYnqZgMoG8q9Ec9k6kGiRpoI4OeJgftIZMvhPgCn\nJGW9IsQcAGaenLzX5US0AhE1ENH3iOhHwLc7durNiWdCkVU77Daz7KEkjGUYzPxHAGcD+BOUJ81/\noFQjO7Ga2wKoxvgwlJ59BpTx93+Tjg4ALoHqYL4iotPFM5mkjOPQfY1DAXxCRLOgRpSHJPl+HcrA\ne22Spw+hVDNI8v2/UB3ddChd+QOeMvgKwLZQo/yXiWg2lNQxE2qezwwAP4Uy9n8J4EwAP02ue9+P\nmecCuAjKM2kGEW3teP+ToUbk4wA8B2V0vtmMT4KIlgOwP4BrmHmq+H0KRfaHud7ZlV+BUVAS2dMA\n/sip+6kZVh5fBWC/5D21He0YKPfmL6Em174gwv8VSp08GkpN9IAR32FQBu33oL7x35FKKFsC+A8R\nfQ1VN0/hGufulIhDu60lRUQDoHSKq0BVlBuZ+WpLuKuhDGBzARwuDJ8l2gBENBTAmswca4QtsZQi\nkTDHAWhyOEOUWMbRnjYMvZzFW4kq4XUiGiEN10S0O1RntVYyOrse8d45JYrBkjYpr0SJEu2EdlNJ\ncdxyFt9OXWfmlwH0oOxM1hL1R9yEnhLfFZTfuoQTS4SXFLmXs7At/9EfWbfAEnUEM5/f3nko0TZI\n7ABFLFVS4juKdicMCi9nYapEbAbAclRUokSJElWAc2xx3a5eUqQWuHsAwB3MbJtlaO7r3R8On/Ci\nps8v7b+hQ4e2ex6WlF9ZFmVZlGXh/+VFuxFG5HIWjyBxD0xm+c7kdHXFEiVKlCjRhmhPlZRezuJt\nItKusmcjmfHJzDcw83BSK5Z+BOWnHlwCoESJEiVK1AftRhgcuZwFM5/UBtn5zmDw4MHtnYUlBmVZ\npCjLIkVZFtWj3SbuFQki4u/Ce5TwgwiYPBno3TsctkQxYAY+/RRYffX2zkmJeoCIwEuL0bvEso0F\nC4AJwml67lxg0SJ3eACYUlqw6o4bbwRuTTbufO01YI01gNZy3ncJfIcIY8EC4E9/yv/cwoXA0UcX\nn58SYVx2GbCaWKO0a1fgEMseaC0twGfJ4t+LF1feL1E7WlrUDwCGDFE/AGhMZmV8/LH6f+WV8hss\nywgSBhGtQ0TPENG7yfnGRHRO/bOWD2++CfzqV/HhZ89WKo733weGDatfvkrEgxn4+9+B6dOz16+4\nAhg0SB2XnVV9sMUWwJ57pudaotAkMneu+t96a+AB51KOJb7riJEw/grlvbQwOR+D+B24vCCim4ho\nChGNcdwfTGrfYb0pkJOopAlj8uSwCP311+pfj5xKtD369LFfnz8/e/7JJ+lxSwtw4onAtdfWL1/L\nIt5+G/iP2JdQtx+tIrzkEuDSS7PXSizZaGkBmgp2a4ohjOVZreME4Nsd14qqMjdDbdbjwyhm3iz5\nXegKJAmib99UB+tCQ/LmWic+f76SOGbMSEdVJeqLrsk+b/Pm+cPNmpUet7QAf/4zcM019cvXsoRp\n01IiIGH61AMw3RbuvRc466y2zVuJMJ57zq6Kv+02YOedi5fIYwhjGhGtqU+IaD/490COBjM/B+Cr\nQLAoC74pUcyxLTIioHWzumFoiaNnTzWaKlF/6LIPVerZs9NjHbbhO2N9a1/cf39KBJIwTAkDqGwz\nJdoGXzl6yDfeAC6+2K6Kv/tuYOTI4vMS0+xOAnADgHWJaBLUbmDHF58VKxjAtkQ0moiGE9H6zoBG\nJV5++cr706dXEos+l9fPO0/ZRErUF7aylxg1Sn03KWFowmgsl8jz4oorKlVHN94IPPqo+xmyDM2k\ntF2Wedvj00+BlVZS7cDs47bYAnjySXU8ezbQqZM9jiIJPqjhYuaPAexERF0ANDDz18UlH8QbAAYw\n81wi2g3AQwDWtgU8++xmAEBzMwAMxvLLD87cP+EE4C9/Ubrv449XXlVA2qjMTuuaa4CddgK22Ua5\nFZaoHptuCtx0E7D55tnrIcIYPBgYNy4rLeoOrOy8/Dj9dGDXXYH11kuvDRkCDBiQGrd1G/BBkk5T\nk/IqLFF/vPcesP76wMyZ6vzEE4HhwxWBAMBDxsp7o0Zlv40kCX18333AqquOxMgaRA8nYRDRGeKU\nxfUkE3x51alGQpITMz9BRH8mopU4uz0mAOCFF5oBKMI4//yshHHuuYosAGWzOPlkpQcH0kK2SR6H\nHgoceWTpRVULmIHRo4GxYysJQ0sLNsLQxMCcHeXqTs6mknr2WWD77UuViYZNYli8WDmFdOkCdO+u\n2oIZnkiV4b/+5ZYwZsxQcfXqVZ+8L8t4911gww2zdf+FF1LXckBpQSR885NaW9U3PeggYPbswZmZ\n7uefn2/3Ap9KagUAXQFsAaWC6ge1WuxxADb3PFcYiGjVZJFCENFWUDPTK8hCQhNAhw7ptQuFqbyp\nKatu0uFNQ7fuxFxiXok46HK0deIuCWPYsPT7EWVtHPp72SSMMVZfu2UXNsL46ivlFHL99ep8xozK\n8Pp/l13cNoytt1YjYBOmh1uJ/NAuzEClfa+5Gdhuu9TmquHzXGttTe9PrtH67JQwmLkZAIjoOQCb\n69F+ssfz8NqSVSCiuwFsD2BlIpoAYCiADkn6NwDYD8DxRNQCtaf3QaE4tfpCd0Kff56939iYdTXT\nI1Zz8p5+vmPH+PcpUQmf2sl175130mMXYZRG7zBshKE90vQ9WY7yWmurGkQdeGB6X6qwxo2zf9PO\nnZWd5Kc/rS3vyzJkfdcDWf1tHnxQuUCvvHIaprExJQTmyu/OnCWMta1K/TjEeOmugqwb7aLkWs1g\nZu98Dma+DsB1eeJ84w31ryvzzjtn759zDrDDDum5bgRPP22mrf6XW06puHbbDdhqqzw5KQG4JYy+\nfVMpwux4ZINpbY2XMEpVlCov7VXmI1Vd5rZytBENkHV/1moOicsTJfXUqXF5LZHFuusqd1hZ300J\nQ59LKWS55dJ2sckmilAkpIQxbVpteYwhjNsAvEJE/4Bycd0HyT7bSyK0ukkX7KRJlWFkI3EZ8aRK\nqrkZ+OAD4M47C8vmMgOXFCFFY1ejANQIS577bBjLMu64A1hnHeD115VTB6DKTtsrzAmqNvdkm9Th\ngyl9a7tgly758l5C4b//VcZrOTA1VeX6u0ktSceOqSpQq2Xl4Gn2bCURArVPuoxZXvwiqH0oZgKY\nAeBwZr64tmTrB11wuoOSPvwaeQkDcI+6SmRHMLZ78h+olATM0ao5woqVMJZl/PznyrlDqmDXW09J\ncvvsU+lw4HL2ePnlOO8poJIwevRQ/6ZLe4l4tLZmv4mpktJtQdpWO3Xyf7NTTgG2/f/tXXmYVMW1\n/53pmWEwKIqyyGJQUSKIPlxxR41KNE9j3LOgRo1JXIgvMW6JGb+4xGeiEZf3jEYTY+K+RkVFkYif\ny3sKLogmKm7gElcEgccs5/1RfbjnVlfdW7dnenq65/6+r7/uvt19b3XdqvrV2Xc0r196CbjuuvLb\nF5JLagMAHwK4E8at9ePisV4J6bikYDC92Pg6Ojd6h+OMM/yLhE0YnZ2lO9gklZQtYcgEyiph68hW\naAAAIABJREFUuDYO9YYNN3RvbBY7ihqLekkTPREwaVL49bRjCRDtenPVYPmwN0hZJQwXtIH8/PON\n52e5CJl29wO4D8C9AB4GsBDAjPIvGSEtl1TxO9OJ6NVi8N7EtHP6dk4aLqO3jdtuM885YaTj2Wfd\n6VRuvdXsbgHg/feN2+tHH5V+zzcp5DPXBMoqYQwcCDz5ZLbf1Bp8ublcC7gQhu77rFK0LWGMHGme\n81iN8uEjjOeeiz4H4humNAlDz5V11ula+0JUUpsz84TiYxMA2wF4Ku13gUjMJUVE+wIYU7zu9wH8\nV9oJL7rIPHd0APfc4/6OPp4mfruSd/3rXyaHSw4Dm5zb2oyXzKWXGl9+wKhLdtvNbVOyd0f2hHGJ\n6EmEIWPAhous6gmNjeGLfiUI48tfNs85YZQP28nDt5nS8S9pEoaeK4MGda19mU2HzDwXwPZdu+zq\nc6XlktofRQN7MQHi2kQ0NOTcnZ3AAQekfy9tcMtiJZPpqquAoUOBXXcNaUXfgAxiGbRPPgncd1+p\nSzNgJA0bdvLBRx6Jn1tPIFGhuFRSspP+2c/ixyTgqd5VJT41net/S5/YKqkssAlDzpUTRvkIJYz2\n9sg43q9fMmHoTa84JJSblDDEhvET9Ti1GDvh0IpWBCMAqJpsWAQTPJiK0A5JkzBswpAdM2AmSL0v\nQj78279FdRGkj846yzxLFT1XkJB2BxQkDXZbJaWN3t/6Vrpt4s47o1oa9Y6GhvBFX9zIyyEMGfO2\nDUMWt1WrjF3roIPCztcX8MILYf2rN0gPPliaeFDXKRF1eXNzuErKTomUFSEShkR8DwDQDGPLCNi7\ndxvsbg5aokNLSqbthpKIZ7PNgOOOC7tOveH554GDDzavpa9F5ZO003T1ZxbCEHViQ4PJyCm6XcAf\n2Vzv0LEuWaWEkIVDVE0CiRK3CaOtzexm29qMC/odd8RL8PZlhPaDljCmTIlSzwu0hCGEkUXCkO+V\nSxghcRgLmPkWfYCIDgFwa3mXzITFAEap9yPhlW5a1evJ6OycHHSBNAlDL1Y33hj/7PXXzeOaa4Iu\nVbeQPpKdZ9JgdBFGUj0Mn1vtWmuZZ5EwTj3VSBM2ktx56wWys+/oyE4YaTaMQYOManf69OiYFFp6\n4QXjovmd7xjyaGszKo/zz49sVbfdBpxySrY29WWsWGH6T2BvulyE4ZIw9FjXEsaSJbMBzMa555qo\n/KwIkTDOcBw7M/ulysI9AKYCABFNAvAZM3vSbLWqx2R0dEQFegT2eyBcJcVsVCBdMRLWK+xo7qyE\nkUXCEIlBBrt4s/3mN+7qia6snfUGGZOu6Os0pO00OzpKPQWJTFQyYFw0JaBVCEM7NuTxMm48+WRE\nBtrW9/LLUbYKoHR9krnQ1hYuYeh78OmnkwG04rvfbcU3v9maud1ewiCirxHRZQBGFF1bLys+/ohu\nqrhXtIc8AWAsEb1DRN8jouOJ6HgAYOb7ASwkotdganL8KPTcnZ1GlD799OiYyygYqpKSRVHvhut1\nAfKByFQy7OyMgrSAbIThcr9NkjBswvjsM3Mf5VhaZcVQ1WQt489/Ns/nnefue3uc6gUkbQPU2VlK\nGMzA+PHR+y++MLVmFi0qjcexCWPmzL4RE9PZGU+Gafftjjuafn3uOZNyXmB7ZdrzRfquvT1yOmhu\nLk1GqOEi7VmzpBRENiRJGO8CeBbAyuKzPO4BsE/2S5WCmY9g5uHM3MzMo5j5Wma+qph4UL5zIjOP\nYeYtix5aQRBVxhCV9crVcaEShixQ9e7Ln4a//930o+Tpb2szcRhAz0gYq1aZmAr5TVJU8Y03Rmky\ndPs0Bg82u7paxg9+YJ5XrjSRvDZs0tTeTWlGbxdhAHG3zrY24zL92GOl9+Okk4Attoje77238TSs\ndzz4YPx/+2A7gdibJ1+5aFslNddaGfVYd93XlSvLk/68hMHMzzPzHwFszMx/YuY/Fh93MHOvNiUS\nRcYjPTlcHZdVwpCFUuDq9KOPrt+64PbuUBNupWwY9kI/cGA00QYPNjtcF+64I/5+/nwjnbz1VuRW\n+tFHplZHvcBVF8Eeiy0t7s+yEIaWMNvaIndpO48Uc7TTlvPXY66pJ56IE3Ooa7Htem4TiG8uaQkj\n7Vquc3Q7YRCRGLXnEtGL1uMF3+96AyQltk0YLpVUqNHbp9pwnfOPfzQiej3CJgy94FSCMFzEu/ba\n0W/eesufStve7T7/vGnj6NHAhAnR8XpawFyZYu0+1wSQZsPo7HSn+Nd929YWkYFP4vvpT6PXmrDq\nBTvtBDz0UPTeXhd8tiWdPh4oJQzfxnP58ug+pi38rjGxYoU7KDkNST+ZVnz+9+ynrS58EoarY0MD\n93z2Cl+wVLmBMb0ddn/pBUf6qhyVlHjZ2NCTULDJJnHDqmTitGF7gWjSX7kyuqcuZ4haRQhhZFFJ\ndXSUus8CfsLwpdL57W+j17rkbj3B55mkQeTOEiG/9Rm5bWijdxphzHAkcqqESurd4vObrkf2S/Uc\nRMLQej6gPBtGORJG0vdrCStXmrKOQOTqZ+82XXWEy5EwfDtTV3nc3XYzJSsFvntgL4D6nhQK0b13\nLYjVxksvuWNI2ttNTYP99zfxEPZGxmVQtseizIkf/tBfB1r/1rUT1WSsCcOWHlzSRJKBtpahx2HS\nYvz226XHRJL44gtTe0f6N8mxxjfubw0IeKiESmoZES31PLrFz4GIphDRK8Xkgqc5Pp9MREuIaF7x\n8fOw84bbMEIJw8f0K1YAl11WerweCGPxYuDmm81r8b13BWsJyiWM//zPbOLxuuvG3/smjt0Oewco\nC5ccf+cd4I03wttRSWy+uYmBsFUS559vHDn+9jfjzhoyzuw+l75ec804YbicDzo73QuLLWHIPbAl\nDNdv65UwCgWzxrzwQvJ4dsU/yDq0bJkh2ZBszLKe2aRy6KHu7+tzVkLCGMDMa3oea2W/VBxEVABw\nOUzywXEAjiCizRxf/TszTyw+znV87ji3mzBcNyHNOH3DDeY5aWKefHL0Wm5eraqkrrkGuOkm81oT\nrHjFhBBGUp/6+iVL/iF71xpKGLaEIbtxOb7NNsBGG4W3oydg95cY64GojGrWc8hi1tgY3zD5nAdc\n/euTMFwuuDZslZQ2/La2+pOG9nbIAvzaa8mLsUuaFpvcF1+YPgwhjCR1uX3Pp08Hjjwyet/thKFB\nRFsR0TQiOomItkr/RRC2A/BaUcXVBuAmuFOOlBUeJ0ZvvcCV00EymEMlBvne3Lkmn06t4bjjgNNK\nZL0ostqVmVZQroQBhKWRP+oo82xLij6DYhphyE5Xji9Zkt6Gnoa9GNj9HbIx8RFGU1OcMJYtK90Z\n+/JT2YZz+c5a1lbSlTtMt+fdd+NxCOecE7kJ1xpkfXFJZboPXRKGrk8SKmHoOt6+8wkaG+Obi4oR\nBhGdDZMxdhCAwQCuI6JfZL9UCVyJBUdY32EAOxZrYdxPRONCTtzUFO2+0iSMUHR2GndOHy69NPoe\nYMpV/vrX5V+vmnB5DUk+IT3ompu7jzBCBq9IIbY04rqvYsPSsFVSImHMm2cmdG90hdZtbm+PG7W7\nShi2hNHS4pbeXP2r79fll0cLYtIccbUntLpfLUAThvSZzy5kQy/woRKGzAPXNWxpsaGhewgjRHP8\nHQBbMPNKACCiCwA8D+BX2S8XQ0ic9FwAo5h5ORF9Dabi36bur7Zi112BjTcGVq6cjI6OySVpDUJu\nQqFgBvRf/mKiIcXo2t5uvGl8u9CZM4Fp06LBoLPa1ir0rqhQADbdNE4Gzc3dY/QGwu6NLC42YbgG\nfmNjspdUY2NEGGKElLYde6yp32En3asGpM0/+YkJGrUDR7tKGPqzzTc3tes1RC9vw75fshiFEIa+\nDy4jfW/NoDB3rlkDNvWsQC6S6Ogoldpcuec0YayxRvp8mDQpmg+u/rIlu+XL5R7NBjAbL79cXmLO\nkD33YgB66rXASANdhZ1YcJR9XmZeyszLi69nAGgiIk8JkFbsvnsrrr22FePGTV5tw9AqqZBFSb4z\neXI8OOmRR9ylLgWyaNWq7ULDtytqbo7vUiQzqf27SkkYX/2qqUMiE0V0wT6VlC2W64WqoaHU6C3n\n+sMfgHvvTW9PT0DafPHFpWRRroQhfa3nhkiL9n3wqaT0XFpvvShgMAthvP++SZMP1EbOr623Bo45\npvS4XQJB1+VesMA86/8kubc09AI/cmR6TrCDD84mYSxbBlxyCQBMBtCK9ddvxc47tyZfxIEQwvgc\nwEtE9MdiHqn5AJYU80pNT/5pIp4BsAkRjSaiZgCHwaQdWQ0iGkpkuo6ItgNAzPyJ74QyiCXXkM3u\nIYQh329oyJYCWIuj9QIZiFttFUX86j6x+6jSEsZ225nUJDJRNt44/LdA/N50dpYavYHoPra3G8mj\n2gkmkxZPu9iOD3b/aAlDUChEhPHaa8m/t49JaVYgLCjPlZet3HTbPQHmaAEeO7b0c9uWoAljyy3N\nc9p90gv8Bhukj+mmpmyEsXQpMEIp/J96qnJG7zthstPOLj7OAnA3otxSZYGZ2wGcCOBBAAsA3MzM\nL+vkgwAOBvAiET0H4HcADk86p3RyoWBUR+3tpZMiDbLramjI5rkj16klwli0KHlBlIE4b55bwpA+\n0v/92WdNtT0ffHYCPUFE5PcZt0XCkMp8oYu6nrSdnaVGbyD6L21tPe9ie+ihpT70nZ3A7Nnu7zMD\nU6emn9d2KNBGb0GhYO5loRARsVwjTcLQbs4hc0z6W3/XpdrsLZLG7bdHwZ2uqHedLRgw7dZjavny\n9PxZ2tV42LB0wtDq4BCVlMuVuRxtSKoNo5hPqiIoqplmWMd04sErAFwRej4tYUiRF93xIQuLDIi+\nIGG4Aog09H8Re5DWOUsfNTVFk2bBAmCPPYwu3NaHA35HAL14yILWv3984Mv9k4kirr6hfa7jDDo7\no12Ynjhy7fb2npcubr3VuM0eckh0rLMT2H139/eZTZK7NNjj2CVhNDZGhGFfo1KEkdTGJUuMSrg3\nkIZOTmn3xQ9/GJUE1kG++j8++GC6ilMv6P37h0kYSTYMmwzEE+2EE0zdmHffLa9vQ7yk/r0YNPdp\ndwfudTfkZupBq29wiOpCE0YWCUPOncTay5YZUfD443tHxk4fIbpiSVatcksYbW1RnzGbQTxkSHaP\nNFeUrC8fzyGHxNPWhxKG9sjRhZl8KqlqwCbxpP8WOuFtTySfhCEV8+zrp6mkhDDuvz86fsghpiKl\nC3YGaKBUtSmpX3rDBkynpyEydrRvf9u8/+//NokHgaitnZ3xe5OUK02gCSPES0pLGGl9tMsu0Xy5\n/HJgzz3N63IkjJBp/TsARwJYtzsD9yqBJFc2wERku1JNaGiVVBYJ44MPzMBIunkXXgjssAPw+99H\nElA1EZJ4TvDJJ24bhhCJfH/VqnC3QA39/bRUH6NHAxdcEB0PHfi2hCGk4MoBpGMLOjsj9VelYQdG\nJo2nLIup7fEGlKpr3323NBdVZ6dJ1JgUVzCo6IYycGD0vX339Vd0+/RT811937Q0+uGHxukESE59\n31N4/PH4+zlz4iSibRfynOYJZkMHM/brly7dumwYRxzh/q4tsXQluDhkWi8C8BIz9wKuT4Z0ygYb\nuD/fdVdTISwJ5RDGoEFmQdl1V/8kfv31eMK83pC/yP5/di0F/V8uucQM5DQJQwgkK2G4pEK7L32T\nqFyVlCvti0sl1dpqPLR6AvYimyRFPP98+Hm1MdpHGIDbs2yzzZJjWtZc0zw3NUXnKRRK78thhwG/\n/KXxNvz88/g5d9klfn+FuFyBf9WEtNG1APtsGCGZq5cujewkoRLG1KmmT+X6V16Z3Ga7vZUijNMA\nzCCiM4joJ8XHf2S/VOWhxeGLLzavs+qhtUoqJPoYAA46yDw/84z/JowZA1x7bfS+NxCGnrALFhg/\nfA1X4Z0kCaO93dRvbm6O6iMk4dFHo9d6sGcljFBi9xGGVtloo7fgV12NOFJ44IG4Os1GS0s8s2tX\nVTLvFENjNRFJP4a4nPsISx/Xld+0HdFu+4gRZoxlWbCqTRi+MegiDJGSbAkjlDAkYDYk0rupycQL\n3XRTdH3fb+zjNsFlQQhh/ArAMpj4iwHFx5rZL1WKtOSDxe9ML37+PBFNTDrftttGr9dbzzyPHm3c\nMUOhCeOKQHO7TocQehOSkpMtXmwK+wjeey/snBrt7WZxSoJeFNNSOACGQN98M3ov0aOy8Igelije\nfh8mTox2pxPVnfUlVfMhtJi9Jox9940IUxOG7JBvvz37ZuNnP0svxvTb3xrVpA8tLfHqaV0hjF12\nicaZy91Vj8GsQVz63sj9T5MwiOJEEmInSiOMt97qOqkcc4x/d758ufkv8h9ddlLpixNOiN53hTBC\nJAxN9nJ9aZOrBrtGpQljfWb+JjP/kpnPkUf2S8URknyQiPYFMIaZNwHwfQBezT+zUQkJhDD69weu\nvz68XVoltfbabr9rGzpYyXUT/va30mNz5hgjoQv77GOITjB8OPA//5PejnffjYILH37YpElOgiYM\nVyGkpNKeQKlKSqLgQweiDp7Ti5cM8K2KWctE9eJawK++GthwQ/81JJEiEBHGttuaSeWSMGTSLVwY\n9j/22CO6NxddFE+TUQ5aWuK/L5cwxowxalIZz5owXBLGhAnAfvuVdy1NGNq13UUYvnriPqSRwejR\n8TK8SfjVr9xz7tpr/TbFL74A1lmnVIrVqTZcWZH1fw9Zf7KqpPT9lPkq98HO59XTKqn7iahbanhb\nCEk+uD9MHisw89MA1iaioSEn18FEWaAJAwhzE9Q3aObM0s/339/9u1NPdR8fPtwdeCPxED5MnBh5\npoSoafR3XIPHV0dBYKukZJcaOhBl4K63ntuG8YMfmO9IbWTXQjxoUHI+oh13jF7L9770pbjRW3vD\n6df2gubq+0cfjReoce0mmZNjU/S5m5vj//OZZ0qvl4YrrjCblKamuIuyDU3SHR0mo6moV5nDU867\nJIyGBreEqI93B2EA4TXZzz4bOO+8+DHZ5Og+f/nlKOnosmVGCtZ2OsD8jylTzGvbm9JWSYVg2bIo\nc4E9BlxIIgy7IJhNPt/4RtTOrAghjB/B2DBWdrNbbUjyQdd3gqhgwgTg1VfN6xDVhuxSbMII+a0m\njGOPDWmdwYIFwN57lx6XHbM9ELfaKnnh+eijaMduL9riqaE9NvSEdU3epEptgOmj6dMjqUbqnact\nBBddZLycZGAPHpwsYQgGOZLC6EJILuhdtEwQIQyXhKFrttv/I4ubtcbHH0dlZH3jSbua6gn+zW/G\nv6eD6lz49FPgRz8CvvIV8z5UJbVsmUntfttt0bGk0rUulVRzc7JKCojvzEMSD2rC6OiIUoksWRKp\n/1x1zEMh6X+IjHvsllsC48ZF0vkXX5hxKv9RBx2K51x3EMbSpfG1x0UYP/95lEnaRRjyG9c81ZCy\nsBWRMIp1MRqYuaWb3WpDw0bsrnP+rrW1dfVjdjE0dswY89mwYem7JSKTK0bYVzrZN6h1TntbBIw1\nNuVfuiQSGTguw/F115UemzHDDFq51sqV8QH73HNGbfbCC3H1mU8lJbAHvf0/OzqM6kYmrPRVGmEs\nW2YMvw0NRprab79St9rXXosb4ZnNfbRRKMSr79lw3fdx4/yEoVNt2xNq773jnmRJucWyQseEJO0u\n0yRenf8MiMaSPv71rxvHEN03rp28Vova8NkwtNHbHvtiw8hCGPo7S5ZE6slp0yJp2lWWNgl3323a\nYkfQz5pl5ggQzY0vvjBzXf6j3Cd9j2zCYM4eFLd0afoa1dBgNgOAmzAEtkONbuvs2bPR2toKoBXz\n5rVmayTCstWCiNYBsAmM4bvYSH4s89XiSE0+6PjOyOKxEphOcGPttdNVNERG/F+0yHj62GkoNA45\nxNggZJIlJVwrh8W1P7eNO+8sPbbvvsCf/xwNnP794ykmRE2iDdGffx4f6LqdPhuGXtRffNEQrAtJ\n/7m1NV4nRER/6e8BA7A663AIGhuTJ6c9eW64wezA3nijVCW1557A0KHAP/5h3tvEN2eOIefx4w1R\n6jbKQiO/6eiIFvcQm4Z28U3SX/sIY4MN3FlQZREaMiQ6NnasKaZz113RZ67MvKNHu113t9vO7SUV\navSWuRgSYyH9smpVFCDY2RnPGp0lXuqJJ6JNoY6gt++RDuBsbIzUseL0oe9RkoTRr18YMS5bBqy/\nfvJ3mKPrlksYkydPxuTJk3HOOWYcz5+fzRwdEul9HIDHADwE4ByY3E+tma7iRmryweL7qcV2TALw\nGTN3QQA1WGcd/2cywXyE8dRTRkes9YRrenzG+vcvT43hioQV2PpJge3lkmaTGDgwnjUzRMJgjnZl\nG27o/29JhDFkSGnZXD2gX3vNlGwNhb2Aiju1wJ484vrpkjDsXbGrT3wGQ0kyJ33iciJIIrauEsag\nQcBee5Uel3NpdZ4Q3ahRZuF59dW4HUbgcyawJTqfDSONMEIWUk0Y9jHf++6Azj5dKEQaBfE89OXB\nAspXSY0dG9leXWMllDDSVFK6nVkRYsOYBmOgfpOZdwcwEUCXa5OFJB9k5vsBLCSi1wBcBWNP6TI+\n+cQsTC7vJXtC2oN6++2Nzl2ThI8wiCSlcBxpbqAuCcOnnxTYhOEaDPYxXe4zxIbBHHc79kHONWGC\neRb9PeBvvyCrh5F9v045Jf7eJgxRm+gCSzLh+/WLq4SSFiJfjI6cS+96pd/32SfSe+s+0dfypeIQ\n+AgjVFX1/PNRDq6ttzaBemut5bZX+AjDluqyeEnZEkaSOhdw34NKEIS9cfERhkC7T7sqO2ZdjD/6\nyKhnJW7GxsYbAyedVB5h+OZUpbykVjLzCnNhamHmVwAEOJumg5lnMPNYZh7DzBcUj11lJSA8sfj5\nlsw813+2bNh4YzNp7c5OIwyBXoh8u/7ly42hyoarpq+GizB8u1PJqBpCGPtYvm4yYJhLB/2iRcYr\nS0MThu4nSb4mEFWV7btu/87Gs8/GVSchSNP72p+LhHH33UatIUnc9twTOPzwOGHY/x9w3weXPttF\nGNpeZTsvaAmjHBvGd77j/w0QSRhZcmQdeaQ/55nuB52bSqukRo4sdWawCWPnnZMLVekxqo/JOcpJ\n0e1DFsLQsCWMDz+MR32HQm9C7N8dcUQ8R5vedPY2wninaMO4C8BMIroHwJvZL9V78dBDRrcOmB2P\njufwqV1kQO2/fzxbZwhCCSPkho4rFq2VQi32OZKgcynZyd8OPthUeQOMjURgE8b228cjoUeONDsh\noJQwFi82i7IPtldUCEKcGTT0ovbAA+b//N//mXve3BwnDNsFE0gPepJ+bGszNo+k79rp1uVYqIRx\n9tkm3ubZZ4Ef/9j/G8A4DHzlK8mGbBsDBwLf/370fpttotd20Sl51kbv++4rDTqVJIeA6ffGxmSC\ndM2FAw6ICDctXuGee8LGFZGxy+l2ynWTCGPnnUvXiF//Gnj66eh9WjoiQZLUYEOTgj2+ZN7Nn2+e\nXX00fbqpKpkVIV5SBzLzp8zcCuAXAK4B8I3sl+q92GuvKNCvsdEU6dEoFEp353IT7r47nQBslCNh\n+KKfxXBo++snkY348svE7eyMS1LM8Z2o7GaYS92OGxvN4JV26fPYhDF8eFiBnSyw05mkwc5zJRKG\nqFM++CBMLZZGGG+9FW08fPfCpZfP4iX1ve8ZyShkQRw2zMQXuFyTQ/G//xvd5wkTTN2ShQvd7SsU\nzOZL5pXAljAKheRF35VRWCeClPbcfrv79/fd55YUbRDFgzxDCcMuU6zPJ+3TCU91RgO7el+ShCGQ\nvtNjxCdhDB9e+l3BSSf5nVaSkClFHDPPZuZ7mLlMj/TaxIABpTprXwr1ELgG36xZpTmUfIvSAQcY\ndUESXO63AknNrAnDHvR65y4LvyYMV4oEIE4YIfaOrmKttUpjFZKg9ezyXhPG448nq22YjYHSF0Uv\n/aj7IY1cgPKM3qHqmP32K29xSMLIkcabzLZz+FLTA5EEIv915cp0CaOjw6hH7U0cYDYfcq9cZU/l\n9+UglDD69XMThmsMHX54lIX3oINM1upQSB+5bE0TJ8Y3YuJ+nlbCuBxUcCr7QUSDiGgmEf2TiB4i\norU933uTiF4o1uMISI5RGYwb53e7E2RJw+wafBJkOH9+lGrCNdiZjZgtBnvfYHis6PScZPzWlcKS\nCENeM5f+b1sl5CKMNHfBrkJyhYXoi30SRmOjuyaHCwsXuotDAXGVFGDsJD41kIswXn89uX6Cblco\nEd97b1it7e5AUoYErbIC/BKGFIUScjn2WPemwPa2cyHU+Gz/Xte3SSKMpqZS+9/227tJRI/P224r\n/d8u1/Zzz423z1VU6uKLo4DcN96IipRJ/9Q8YQA4HcBMZt4UwCPF9y4wgMnMPJGZM6QQ7D5cfLHR\n9/kGlCBL9lnX4BMvlhNOiOIAOjtL6ysLRE2UNhiS3GW1hGF/rxzCGDHCrZK66KJ4avfuhk1ISWqX\nhoZSwli1yjzrvvTdT1cfaMhCIc9LEvwJ58wx40sXc1qxIp7V2NV+QXcafLsL0iaX56BNGC6iBiJy\na2mJx7PY0FK/vnezZkUR0eVW7NNxGIWC37PRpZLyVcPzjSlRq2rCkJips84yz+LZ5kJDQ3Tu0aNN\n2MA770TBibareVdQLcJYnSOq+JxkE+lGfsyOU04xxr5Jk+KDxh7EWdQuLsKQCaYX1s5OYJNN3OcI\nJYyOjniGWcB4cWjYhMEc9/xykYfr/VNPxVUHjY1Gv96/f2WljO9+Nx51bUc7a6xaFb9Xzc2mj2xV\nlc+Y7sqx5PKSEokzKdbgrLOMY0FTU7STzIJKqvqywK7cB5gASBsuCcOlkpJzCGH4/qfPRfvoo6NY\nnnIljEcfNVHeHR2mjUkSRihhDBvmJrBJk8yznoM6hc9774UnWBRIPMduu7nvRbnwDjlBCEqZAAAY\nwUlEQVQiWqZyR9mPruaSGqoC8D4A4PtLDOBhInqmGEBYNZx2Wjx9Qsju7je/cR936SFlIEnkM5Ac\nfS0DOG3REJE+CS4JQ+fn0RLG4MHAf6hqKLofRo6Me5g1NBgPnkrXxiaKDHyAyZXjmuDvvGPUBbaE\nIc+hu/cQCUPUSq5cYQIZT8xmgeotBJAVO+8MPPmkeS3jSHsO6gwDLpWUjzAko3BWCUNLdV1JD/+P\nf5Rn9PYRhh6jGvI/9Hy/444oy+2wYb1HmvQO0WIOqTU9j9RcUkUbxYuORyx3KzMz/HmldmLmiQC+\nBuAEItoly5+rJFw30E6R7NsBucRbbQgU2PmggGjn0RWVlA2bMDo64m66mjAaG+MFfpLcWrN6j3UX\nzj/fpBm3MXJk6S5XE4bLD//AA+PnSMsTJJLOpZeaZ1u609AbkI8/Tg9qtFGuuqW70dAQ7ZJlPOg+\n1uThkjBsotTBaaEShmsevPxystFb/yYpViErYTQ2lh475RSzkXGtGy0txqZx9NHRsQMPNJJzb0No\nLqktAewKs7DPYebU4pDM7EhUsPp8HxDRMGZ+n4jWB+BMH8bM7xWfPySiO2Eizue4vqtzSUm+lErC\nNYjtxdMXDewafGmpOX760/hn2tXV/q5GiJeITRhvvWVSir/6qtkpa8Kw4dv5zJuXHJBVaSQtpi5b\nhU8ldcstpcVqQvr0qafSv6NT2H/+uVFLZnGe8AWMVhNjxvhTkttkvWKFUR/ac0mrpJIi320J48IL\nTdEowbhxUWZWF4YMiSTpeL4lo/uXe10o+HObFQpuCUOyRss4FDuCy/mgpSVKK19pzJ49e3Vy1nKQ\nShhENA3AcQDugLEn3EBEVzPz9LKvanJEHQngwuLzXY7rrgGgwMxLiehLAPaGyWXlRFLywUrANYhd\neYsAYKed4vl3bJWUSyUEJC9Mdrpl33dDErO5FsF99jG7Y5nUPvgkDElDXS3YBDp1qvszH2EIbEIM\nJYwQ2DVPskgY06ZVT4JLg89AbBdQWrrUJE30qaTWWAN4+21/f9uEcfrpJshUbxaSVFJ77WWSUfrO\nK84IhYIhnoaGUgLSgYgCl11D4HLI6O7YpCTYm+lzzunm5IMAjgWwPTOfzcy/ADAJhkC6gl8D2IuI\n/glgj+J7ENFwIpKkCcMAzCGi5wA8DeBeZn6oi9ftNrh21q40FECpx4g9QNrbw5L/ua4vk8mnekry\n0tHXsX8/alS0mCZ5nPQW3aqN/faL21oOPjh6rRcgnWlVE4b0vataWVfrbPuQhTCqXeu6HNgSxvvv\nl7o56/f9+xvV3l0l28nou/rcgKmopw3GSfdK4pGAeEVLnXhUCIPI7ankkzBcNgzAnfi0Jwmjqwis\nqYVOz+uywMyfAPiq4/i7APYrvl4IoMr7VD/SCOPRRyNjn00YNrGUSxh61+zbhYXUE3YRxsiR0YSU\nAe0ijNCqbD2N9deP21o0XIRhe+skSRH2ZzapnHii24aShiyu2fVAGIsWlRL1qlVxCSMJPqO3Rsim\ny4a0Z+XKuNHd9f2sEsamm8bf779/aSLK3owQCeM6AE8TUSsRnQPgKQAJnuJ9A0mEseeeRg8qA9rO\nyGn/dskSdyGlpEWrUDA1MAQ+CcNWe7jgIox11y3d8daShGFD7+J0v8q98UkYNkJUUnvsUV4b+xph\niKFf8qF9+csmyDCUMPT49Nk5kmxJvt+4JAzAP+dt12nXMcE228Tv3d13R9URawGJhEFEDTDqoKMB\nfArgYwBHMbMjaXffglTz05AJLyH/PpWUvSu/5BK3LjXJw8kevL4FLqQeh4sw+vcPW8B8hv3ehLlz\nga8qeVYv+GKEtAkji4Rhw5fuPg1ZVFIhG4HeCHuRbm6Oci1tv328+mJ3SBh28sOkttjnsiUMn43L\nJgetknIhrbxBb0aiQoGZO4noCmb+NwDP9lCbagJnnx1lcxXILlYGmI8wCgVDEJKS2uUZ4xJ1NewF\n3lWwBggrUnP77aWLoIswbAnjiSdqY3ekE74BfsLQ/ZckYdifpZWyDYWdot6HqVOT4zt6M1wZErQd\nCQiXMLTUWE6sTxphvPpqPLuxS/3ar1/p/E0jjFpGiErqYSI6mKjS4Ve1hcbGUhc5IQYZWDIRbPdH\nuwCNEMNhh0XSicuXW8MmE58tIUTCmDatlIBaWkzCOk12NmHssENy9cLeCk0Ysku14zB8hKEr9enz\naYkvxN3Vpd5w1XAXbLBB9PpPf4obbGsF66xTurBrqcpOVplGGK4A2CxII4wzz4zfW9f1WlpKpb2+\nThg/AHALgFXdFelNRIcQ0UtE1EFE3uTMRDSFiF4holeJ6LSuXLMnIAuFLWHYImihEF90ZPEfNswY\nwQAz6A491H8tW4ftsyWEDlyXSuq66+IR370lUKyr0H0vi0aoSkrnfRJ0dprFTVyJQ9R0WT1jfNlx\nawULF5r0Fkn1p+1aK2mqG1nAv/3tdM81TVRSCiAksj6EMOy5mLbZq2WE1MMYwMwNzNyUJdI7BS8C\nOBCmVrgTRFQAcDmAKQDGATiCiDbr4nUrCtmN69QGQOnAb24uDZQD4lGvad5Htrus7/t2nQwfXBJG\nY2Nt61t9+NrXIn966e8110yuMSBwEYYck0XLrgHhQlbbT62T9YYblm6UALeEkZUwxo4Nizeyr+Mj\nDHGfHTw4ThguyVGCCzWSvKRqHamEQUSPhBzLAmZ+hZk9CaJXYzsArzHzm8zcBuAmAAd05bqVhgwo\nWyVl7yb79YsP8Icfjn4ng/Ozz5Kv9bkl4/kkjFD3Thdh2Kj1RUvwpS8BN95oUonLomGrF32LSXt7\nOmEMHJgeuCj9e+ut5jkpGylQP33vI4xddjF1XoBwlZR83tKSTQWURhjf+56pxjhkSLqEIcSvtQFr\nrmnaM3Rocl2aWkRS8sH+RLQugMHF+hXyGA1gRA+0bQQAXRJ9UQ9dt2zYEoZPJaULv2gUCuFuqiES\nRhavm7/+Nf6+ngkDMNLERhtF7wcOdFc4tNHWVrqjFLLQWUHT1B2y0GyxhXlOsnt8+mn99L2vLMBj\njxl3dI208SsL+FprlWZgToJdNdJGY6Nx8W1rixOGSyqUeaITC0pRqb/+FTjqqPB21QKSFB/HA5gG\nYDjiHlJLYVRFiSCimTDR2jbOZOa/BbQt0xTp6VxSLsiAksXGRxi2SkqgJYwkTJ0arzB2113ANxwJ\n4keOLC2hGQpXO3z5dGoZklW2qcldp9qGpLzWEAnjllui86URhiw0MmZks/H66/F+HjPGpGapF5eT\nMWNM8a/2dpNgz0UKocGgQhhrr2089kKRJmFIfYn2dpN2/mc/M8dd90Duo5zzySeNcX/QoPLjcSqJ\niuWSYubfAfgdEZ1cTt6opOSDgVgMYJR6PwpGynCip3NJ+bBwYTT5fTYMV5UuwJ2504WxY+OLls9T\nacSI8gnDxvLltZXCIBS+Osq+BXrZslLCWLHCSBiyUADpnlJTppjfyfVtdaZAdqsXXBC5YdcyiExk\n89tvm/cuwlhrLZNn6oorks8lhJG1mmAaYRCZ+yBz9PmEVKsyJyS784oVJnNvSIaFaqDiuaSYeToR\n7UhE3yKiqfLI3FI/fHunZwBsQkSjiagZwGEwSQt7NXSNY58NQ3YvNkJVUqF+/r6UAxMmJP9u221L\nj/XvXz+7XI3jjgO+9S3zWquafIvJvHmRJ5vgqqvMs+6fm29Ovu7pp5tymkIYsvj5ijOtu2681kit\nQzZVPsP2gAHuRXfKlOi1zKushJGmkhIJQ+IrJPMsAJx8cvy7cv/EmzCpvG49IMTofQOA3wDYGcC2\n6lE2iOhAInoHJpHhfUQ0o3h8dfJBZm4HcCKABwEsAHAzM7/clev2NHwqqaamaJe6leVUnEQYEyea\nqOV99okf97mA+kqV+qr4CXQitnrHBhtE6r0QG4avlreNIUOSdfC2Xny33cxrIQy5fj2SNGCksX/+\nM56O3IYuvCTQ80P6Nymbsguumh32NZqaoutrhwR7syUutStXGvXVzjtna0utIURbuDWAccVCR90C\nZr4TwJ2O46uTDxbfzwAwo7uu29OwCePee82uX2IsPvwwnvY8qboYYMgCKC3Ko6WVBQui3Dy+xcY3\nUY49FrjmGv/16x0ywl1V4MpBkj++TRhnnmlKtn7yiXnfv79ZjGq1Cl8I0jYu0hcaWgKT+ZU1sl76\n1Ne3668fD6698cboM9vhQaScjz+O6m/XM0KG43wAFazIXL+QRUFEYHG1bGoyWSsvvTQ+AEONm3bK\nDk0YY8fGP5PzaVF5yy3d581SvKceIQvEW291z0KdlItLJ7Rjju6T/EbsUvVMGGm47DJgjlUuzUUY\nTU2llRGTkEQY48cbiULfO60hsLfNe+5pJP4DerXDf/chZDgOBrCAiB4ior8VH73eltCbIAPTdrcF\n4oSx9dZh53OlRxfdrp4ERNFCJDuhSy4xmXRd2GGHsOvXK+RejBgRLmEcfLDx2XfBzo2k4ZMk5d6K\nOvHnPw9rRz1i1Cij4tF95Sqv29CQLadZmoShz31cSuUfInP/r7wy/Pq1jBDCaAXwDQDnwdgyflt8\n5AjAxIlx/3wgvnuRRer114Edd/Sfx1VWVNDR4U8j4qtmBgBXXx3/rN71r2kIsWGIgVxw663+7LRy\nn3RNa0EaYYhLptTK7stob4/UrK5UImKkDoUmjAcfdHufyX3QWY4B49Iucy00WWQ9IcRLajaAVwCs\nBWBNAAuY+e8VblfdYO7cSKQV47R2gxURVweRuaAnhC1hbL+9e7ekJQzB+PHRsWOPjX/WG2tE9yS0\n5Lfeem7SuP760mPdSRhy/NxzK1fZrxYhUd0NDZH6RxNGliBVTRh77x2RkYYvWK9/f6MqO+ec7Mb2\nekCIl9ShMDUxDgFwKID/IaJDunLRDMkH3ySiF4hoHhHVvO9Oe7uxE+io0NBFQROGvZtaf33/jtjO\nj5QUTLTRRsC//hXWnnrEbrtFubfmzCl1LgDcBvHx441Dgw25TwceWOqxllTt7fe/N2629eohVQ5c\naTmEJIjKlzCAuF0iKTmiYMgQU96gLyJEJfVzANsy81RmngrjUvuLLl43NflgEQxgMjNPZObtunjN\nquK73zXZaO0dSyhhjB8fvXYNYp+EkYZhVix+Wk6jegZRZEcaOjRKKX7FFUBrK3DKKea9fQ8bGkzh\nHxtyn847z6Qu17ajpHtz3HE5WdgQaVj3S4hKyhWXZrssJ83BLETUFxBCGARAZ2r5GP5guyAEJh/U\n1695XH+9OxdNKGHMmhW9du1OQyQMF957Lx5smKMUP/oR8MtfAhdfbN6HRrxr1WFTk6nznlQBLocf\nLjuDNnoLeWgXWAA46ST/OWVuZJUw+jJCCOMBAA8S0VFEdDSA+9FzsREMU8DpGSJK8VeoTfhKbdrR\nxFokd5GAljDk9WGHhUkeJ58MHH98eltzGIQShmuxGTwYOOaY7m1PX4Odq01LGDvtFH3vyivD7pXe\ntOWEkYzUwD1mPpWIDgIgt+KqYuBdIroh+SAA7MTM7xHRYAAziegVZp7j+mJvSD5YDhYvdhPAGWcA\n9wQ4L7sigiVgzJfOwz724x+HtzeHiWNJqo4ncC02hULfDo7sDsii7jJ665iJzz8Pi2NJkvJrsaJk\nEiqWfJCINgEwlJkfZ+bbAdxePL4zEW3MzK8nnbgbkg+Cmd8rPn9IRHfC1MhIJYxawqGHmgyeXYWe\nGPfeG0UYuwgjNIV6DjfuvNM4Gtgp5m2cfHJp4FmOriHNhqEJY/HiMMLQ39HeVu+/X+oSX+voavLB\nJAnjdwDOcBz/vPjZv2e6kh9OLTsRrQGgwMxLiehLAPYGkO3f1QBcSepmzAC2284YWS+5JPn3Lglj\nr71KP9fYdlvg8ceztzWHQf/+xqUyjTC+/e3arL3dm7H55lEqc5cNQ+xGL75o0vuHEMZPfmK8BwcN\niicyrDey6A6QL0UUET3DzNt4PpvPzJuXfVGiAwFMB7AegCUA5jHz14hoOICrmXk/ItoIwB3FnzQC\n+AszX+A5X3emuupVsFMt6+OAmRAdHSZR2tlnl0acDhhg7CR12j1Vw4YbGrfbvF97Du3tRjo+8EDg\n7ruj3GsdHUaa2GAD853ly6PYGOZS0pBULG++aQol9WUQEZg52LEoScJICkvpUmWEkOSDzLwQQEqh\ny76BkPrP663nTk+Qu2dWBn05x1O1INKDSAENDcDTT5vnUaOMLYIoHkiZNP5DCzXliJA07J8hou/b\nB4veSs86vp+jAmhpMZPBhzRCyAmjMth9d3eEcI7K47LLojTz26norCxjff58kzMsRzYkqaSGwUgB\nqxARxNYA+gE4UAzSvQH1rJJ67z2jn7XTS8jkKBTcxZgE995rqsQdfnjl2pgjR2+GTSR1ulSUhawq\nKS9hFE9GAHYHsDlMTMRLzDzL+4MqoZ4JwweZBIcfXhqslCNHjgg5YfjRrYRRKRDRRQC+DiO9vA7g\naGYu8TkhoikwHlkFANcw84We8/U5wviv/zKeOkccUe2W5MjRu5EThh9ZCaNapruHAIxn5i0B/BMO\n910iKgC4HMAUAOMAHEFEm/VoK3sxfvhDN1l0JSin3pD3RYS+3Bfz55ukmrNmmee+3BddRVUIg5ln\nMrPEVz4NYKTja9sBeI2Z32TmNgA3Aegjda3KRz4ZIuR9EaEv94VU0dt9d/Pcl/uiq+gNzoHfg8lP\nZWMEgHfU+0XFYzly5MiRowqomCdySC4pIjoLwCpm/qvje7mmMUeOHDl6Eapi9AYAIjoKwHEA9mTm\nlY7PJwFoZeYpxfdnAOh0Gb6JKCeXHDly5CgD3RXpXTEUvZ9OBbCbiyyKeAbAJkQ0GsC7AA4D4PQJ\nyvKHc+TIkSNHeaiWDeMyAANgUpbPI6IrAYCIhhPRfQDAzO0ATgTwIIAFAG5m5per1N4cOXLk6POo\nmkoqR44cOXLUFnqDl1TZIKIpRPQKEb1KRKdVuz3VAhGNIqJHieglIppPRCdXu03VBhEVitJraLGu\nugQRrU1EtxHRy0S0oGgb7JMgojOKc+RFIvorEQWk9awPENG1RPQBEb2ojg0ioplE9E8ieoiIkhLO\nAqhhwsgD+2JoA3AKM48HMAnACX24LwTTYFSZfV2EvhTA/cy8GYAtAPRJtW7RFnocgK2YeQJM9oi+\nlGHtOpi1UuN0ADOZeVMAjxTfJ6JmCQN5YN9qMPP7zPxc8fUymEVheHVbVT0Q0UgA+wK4Bp4CXX0B\nRDQQwC7MfC1g7IKuFDx9BJ/DbKzWIKJGAGsAWFzdJvUciqWtP7UO7w/gT8XXfwLwjbTz1DJh5IF9\nDhR3UhNhIuj7Ki6B8cJLqNbcJ7AhgA+J6DoimktEVxcrWfY5MPMnAH4L4G0Yr8vPmPnh6raq6hjK\nzB8UX38AILXGYC0TRl9XNZSAiAYAuA3AtKKk0edARF8H8C9mnoc+LF0U0QhgKwBXMvNWAL5AgNqh\nHkFEGwP4MYDRMNL3ACLKC+gWUczemrqm1jJhLAagSwuNgpEy+iSIqAnA7QBuYOa7qt2eKmJHAPsT\n0RsAbgSwBxFdX+U2VQuLACxi5v8tvr8NhkD6IrYB8AQzf1x02b8DZqz0ZXxQrHsEIlofwL/SflDL\nhLE6sI+ImmEC++6pcpuqgmLdkj8AWMDMv6t2e6oJZj6TmUcx84YwRs1ZzDy12u2qBpj5fQDvENGm\nxUNfBfBSFZtUTbwCYBIR9S/Ol6/COEX0ZdwD4Mji6yMBpG40a7aqLTO3E5EE9hUA/KEPB/btBOA7\nAF4gonnFY2cw8wNVbFNvQV9XXZ4E4C/FTdXrAI6ucnuqAmZ+vihpPgNj25oL4PfVbVXPgYhuBLAb\ngPWI6B0AZwP4NYBbiOgYAG8CODT1PHngXo4cOXLkCEEtq6Ry5MiRI0cPIieMHDly5MgRhJwwcuTI\nkSNHEHLCyJEjR44cQcgJI0eOHDlyBCEnjBw5cuTIEYScMHLkyAAiWreYNn0eEb1HRIuKr5cS0eXV\nbl+OHJVEHoeRI0eZIKJfAljKzBdXuy05cvQEcgkjR46ugQCAiCZLsSYiaiWiPxHRY0T0JhF9k4h+\nQ0QvENGMYnptENHWRDSbiJ4hogckr0+OHL0VOWHkyFEZbAhgd5iaAzfAFKrZAsAKAPsVk0VeBuAg\nZt4GpsDNedVqbI4cIajZXFI5cvRiMIAZzNxBRPMBNDDzg8XPXoRJsb0pgPEAHja58FCAqdOQI0ev\nRU4YOXJUBqsAgJk7iahNHe+EmXcE4CVm7usptnPUEHKVVI4c3Y+Qwk3/ADCYiCYBpp4JEY2rbLNy\n5OgacsLIkaNrYPXseg2UpljnYh36gwFcSETPAZgHYIdKNjRHjq4id6vNkSNHjhxByCWMHDly5MgR\nhJwwcuTIkSNHEHLCyJEjR44cQcgJI0eOHDlyBCEnjBw5cuTIEYScMHLkyJEjRxBywsiRI0eOHEHI\nCSNHjhw5cgTh/wHm8vUk/FOvlwAAAABJRU5ErkJggg==\n", "text": [ "" ] } ], "prompt_number": 10 }, { "cell_type": "heading", "level": 3, "metadata": {}, "source": [ "Versions" ] }, { "cell_type": "code", "collapsed": false, "input": [ "from qutip.ipynbtools import version_table\n", "\n", "version_table()" ], "language": "python", "metadata": {}, "outputs": [ { "html": [ "
SoftwareVersion
Numpy1.9.1
Cython0.21.2
IPython2.3.1
SciPy0.14.1
Python3.4.0 (default, Apr 11 2014, 13:05:11) \n", "[GCC 4.8.2]
matplotlib1.4.2
OSposix [linux]
QuTiP3.1.0
Tue Jan 13 13:31:16 2015 JST
" ], "metadata": {}, "output_type": "pyout", "prompt_number": 11, "text": [ "" ] } ], "prompt_number": 11 }, { "cell_type": "heading", "level": 3, "metadata": {}, "source": [ "References" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "[1] Machnes et.al., DYNAMO - Dynamic Framework for Quantum Optimal Control. arXiv.1011.4874" ] } ], "metadata": {} } ] }