{ "cells": [ { "attachments": {}, "cell_type": "markdown", "metadata": {}, "source": [ "# Quantum Enhanced Auxiliary Field Quantum Monte Carlo\n", "\n", "This work was done in collaboration with the Next Generation Computing team at [BASF](https://www.basf.com/global/en.html).\n", "\n", "In this tutorial we implement a quantum-classical hybrid workflow for computing the ground state energies of a strongly interacting molecular system. The algorithm consists of two parts:\n", "\n", "\n", "1. A variational quantum eigensolver that uses the quantum-number-preserving ansatz proposed by [Anselmetti et al. (2021)](https://doi.org/10.1088/1367-2630/ac2cb3) to generate a quantum trial wave function $|\\Psi_T\\rangle$ using CUDA Quantum.\n", "\n", "2. An Auxiliary-Field Quantum Monte Carlo simulation that realizes a classical imaginary time evolution and collects the ground state energy estimates.\n" ] }, { "cell_type": "code", "execution_count": 1, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "\u001b[33mWARNING: Running pip as the 'root' user can result in broken permissions and conflicting behaviour with the system package manager, possibly rendering your system unusable.It is recommended to use a virtual environment instead: https://pip.pypa.io/warnings/venv. Use the --root-user-action option if you know what you are doing and want to suppress this warning.\u001b[0m\u001b[33m\n", "\u001b[0m" ] } ], "source": [ "# Package installs\n", "!pip install pyscf==2.6.2 openfermion==1.6.1 ipie==0.7.1 -q" ] }, { "cell_type": "code", "execution_count": 2, "metadata": {}, "outputs": [ { "name": "stderr", "output_type": "stream", "text": [ "/usr/local/lib/python3.10/dist-packages/qutip/__init__.py:66: UserWarning: The new version of Cython, (>= 3.0.0) is not supported.\n", " warnings.warn(\n" ] } ], "source": [ "# Relevant imports\n", "\n", "import cudaq\n", "\n", "import numpy as np\n", "import matplotlib.pyplot as plt\n", "\n", "from pyscf import gto, scf, ao2mo, mcscf\n", "\n", "from afqmc_src.vqe_cudaq_qnp import VQE, get_cudaq_hamiltonian\n", "from afqmc_src.utils_ipie import get_coeff_wf, gen_ipie_input_from_pyscf_chk\n", "\n", "from ipie.hamiltonians.generic import Generic as HamGeneric\n", "from ipie.qmc.afqmc import AFQMC\n", "from ipie.systems.generic import Generic\n", "from ipie.trial_wavefunction.particle_hole import ParticleHole\n", "from ipie.analysis.extraction import extract_observable\n", "\n", "from ipie.config import config\n", "\n", "cudaq.set_target(\"nvidia\")\n", "\n", "# Ipie has recently added GPU support however this remains a bit tricky to use as it requires manual installation of several packages.\n", "# Once this is streamlined, we can set the GPU option to True in the tutorial.\n", "config.update_option(\"use_gpu\", False)" ] }, { "attachments": {}, "cell_type": "markdown", "metadata": {}, "source": [ "We start by defining the structure of the molecule, the basis set, and its spin. We build the molecule object with PySCF and run a preliminary Hartree-Fock computation. Here we choose as an example a [chelating agent](https://doi.org/10.1021/acs.jctc.3c01375) representing a relevant class of substances industrially produced at large scales. Their use ranges, among the others, from water softeners in cleaning applications, modulators of redox behaviour in oxidative bleaching, scale suppressants, soil remediation and ligands for catalysts. In particular we focus here in a Fe(III)-NTA complex whose structure is given in the file imported below.\n" ] }, { "cell_type": "code", "execution_count": 3, "metadata": {}, "outputs": [], "source": [ "# Define the molecular structure and the basis set for the Fenta molecule.\n", "\n", "atom = \"afqmc_src/geo_fenta.xyz\"\n", "basis = \"cc-pVTZ\"\n", "spin = 1\n", "num_active_orbitals = 5\n", "num_active_electrons = 5\n", "\n", "# You can swap to O3 which is a smaller system and takes less computational resources and time to run.\n", "\n", "atom = \"afqmc_src/geo_o3.xyz\"\n", "basis = \"cc-pVTZ\"\n", "spin = 0\n", "num_active_orbitals = 9\n", "num_active_electrons = 12" ] }, { "cell_type": "code", "execution_count": 4, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "-224.34048064812222" ] }, "execution_count": 4, "metadata": {}, "output_type": "execute_result" } ], "source": [ "# PYSCF helps to build the molecule and run Hartree-Fock.\n", "\n", "# Define the molecule.\n", "molecule = gto.M(atom=atom, spin=spin, basis=basis, verbose=0)\n", "\n", "# Restriced open shell HF.\n", "hartee_fock = scf.ROHF(molecule)\n", "hartee_fock.chkfile = \"afqmc_src/output.chk\"\n", "\n", "# Run Hartree-Fock.\n", "hartee_fock.kernel()" ] }, { "attachments": {}, "cell_type": "markdown", "metadata": {}, "source": [ "### Hamiltonian preparation for VQE\n", "\n", "Since this molecule contains of around 600 orbitals (which would correspond to 1200 qubits) and 143 total electrons, it is impossible to perform a full VQE with full statevector simulation. Therefore, we need to identify an active space with fewer orbitals and electrons that contribute to the strongly interacting part of the whole molecule. We then run a post Hartree-Fock computation with the PySCF's built-in CASCI method in order to obtain the one-body ($t_{pq}$) and two-body \n", "($V_{prqs}$) integrals that define the molecular Hamiltonian in the active space:\n", "\n", "$$ H= \\sum_{pq}t_{pq}\\hat{a}_{p}^\\dagger \\hat {a}_{q}+\\sum_{pqrs} V_{prqs}\\hat a_{p}^\\dagger \\hat a_{q}^\\dagger \\hat a_{s}\\hat a_{r} \\tag{1}$$\n" ] }, { "cell_type": "code", "execution_count": 5, "metadata": {}, "outputs": [], "source": [ "from openfermion.transforms import jordan_wigner\n", "from openfermion import generate_hamiltonian\n", "\n", "# Run a CASCI simulation for computing the Hamiltonian in the active space.\n", "casci = mcscf.CASCI(hartee_fock, num_active_orbitals, num_active_electrons)\n", "casci.fix_spin_(ss=(molecule.spin / 2 * (molecule.spin / 2 + 1)))\n", "\n", "# Executes the kernel to compute the hamiltonian in the active space.\n", "casci.kernel()\n", "\n", "# Compute the one-body (h1) and two-body integrals (tbi) as shown in equation 1.\n", "h1, energy_core = casci.get_h1eff()\n", "\n", "h2 = casci.get_h2eff()\n", "h2_no_symmetry = ao2mo.restore('1', h2, num_active_orbitals)\n", "\n", "# V_pqrs terms in H.\n", "tbi = np.asarray(h2_no_symmetry.transpose(0, 2, 3, 1), order='C')\n", "\n", "# Compute the hamiltonian and convert it to a CUDA-Q operator.\n", "mol_ham = generate_hamiltonian(h1, tbi, energy_core.item())\n", "jw_hamiltonian = jordan_wigner(mol_ham)\n", "hamiltonian, constant_term = get_cudaq_hamiltonian(jw_hamiltonian)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Run VQE with CUDA-Q\n" ] }, { "attachments": {}, "cell_type": "markdown", "metadata": {}, "source": [ "We can now execute the VQE algorithm using the quantum number preserving ansatz. At the end of the VQE, we store the final statevector that will be used in the classical AFQMC computation as an initial guess.\n" ] }, { "cell_type": "code", "execution_count": 6, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "# Using cudaq optimizer\n", "# Num Params: 16\n", "# Qubits: 18\n", "# N_layers: 1\n", "# Energy after the VQE: -224.3881035525103\n" ] } ], "source": [ "# Define some options for the VQE.\n", "options = {\n", " 'n_vqe_layers': 1,\n", " 'maxiter': 100,\n", " 'energy_core': constant_term,\n", " 'return_final_state_vec': True\n", "}\n", "\n", "n_qubits = 2 * num_active_orbitals\n", "\n", "vqe = VQE(n_qubits=n_qubits,\n", " num_active_electrons=num_active_electrons,\n", " spin=spin,\n", " options=options)\n", "\n", "results = vqe.execute(hamiltonian)\n", "\n", "# Best energy from VQE.\n", "optimized_energy = results['energy_optimized']\n", "\n", "# Final state vector.\n", "final_state_vector = results[\"state_vec\"]\n", "\n", "# Energies during the optimization loop.\n", "vqe_energies = results[\"callback_energies\"]" ] }, { "attachments": {}, "cell_type": "markdown", "metadata": {}, "source": [ "### Auxiliary Field Quantum Monte Carlo (AFQMC)\n" ] }, { "attachments": {}, "cell_type": "markdown", "metadata": {}, "source": [ "AFQMC is a numerical method for computing relevant properties of strongly interacting molecules. AFQMC is a type of Quantum Monte Carlo method that combines the use of random walks with an auxiliary field to simulate the imaginary-time evolution of a quantum system and drive it to the lowest energy state. This method can provide accurate results for ground-state properties of a wide range of physical systems, including atoms, molecules, and solids. Here we summarize the main features of AFQMC while a detailed introduction to can be found [here](https://www.cond-mat.de/events/correl13/manuscripts/zhang.pdf).\n" ] }, { "attachments": {}, "cell_type": "markdown", "metadata": {}, "source": [ "We consider the electronic Hamiltonian in the second quantization\n", "\\begin{equation}\n", "H = {H}_1 + {H}_2 \n", "=\\sum_{pq} h_{pq} {a}_{p}^{\\dagger} {a}_{q} + \\frac{1}{2} \\sum_{pqrs} v_{pqrs}{a}_{p}^{\\dagger} {a}_r {a}^{\\dagger}_{q} {a}_s \\tag{2}\n", "\\end{equation}\n", "where ${a}_{p}^{\\dagger}$ and ${a}_{q}$ are fermionic creation and annihilation operators of orbitals $p$ and $q$, respectively. The terms $h_{pq} $ and \n", "$v_{pqrs}$ are the matrix elements of the one-body, $H_1$, and two-body, $H_2$, interactions of $H$, respectively. Here, we omit the spin indices for simplicity.\n", "\n", "AFQMC realizes an imaginary time propagation of an initial state (chosen as a Slater determinant) $\\ket{\\Psi_{I}}$ towards the ground state $\\ket{\\Psi_0}$ of a given hamiltonian, $H$, with\n", "\\begin{equation}\n", "\\ket{\\Psi_0} \\sim\\lim_{n \\to \\infty} \\left[ e^{-\\Delta\\tau H } \\right]^{n} \\ket{\\Psi_{I}}\n", "\\tag{3}\n", "\\end{equation} \n", "where $\\Delta\\tau$ is the imaginary time step.\n", "\n", "AFQMC relies on decomposing the two-body interactions $H_2$ in terms of sum of squares of one-body operators ${v}_\\gamma$ such that the Hamiltonian ${H}$ becomes\n", "\\begin{equation}\n", "H = v_0 - \\frac{1}{2}\\sum_{\\gamma=1}^{N_\\gamma} {v}_\\gamma^2\n", "\\tag{4}\n", "\\end{equation}\n", "with ${v}_0 = {H}_1 $ and $\n", "{v}_\\gamma = i \\sum_{pq} L^{\\gamma}_{pq} {a}_{p}^{\\dagger}{a}_{q}.\n", "$\n", "The $N_\\gamma$ matrices $L^{\\gamma}_{pq}$ are called Cholesky vectors as they are obtained via a Cholesky decomposition of the two-body matrix elements \n", "$v_{pqrs}$ via $v_{pqrs} = \\sum_{\\gamma=1}^{N_\\gamma} L^{\\gamma}_{pr} L^{\\gamma}_{qs}$.\n", "\n", "The imaginary time propagation evolves an ensemble of walkers $\\{\\phi^{(n)}\\}$ (that are Slater determinants) and allows one to access observables of the system. For example, the local energy\n", "\\begin{equation}\n", "\\mathcal{E}_{\\text{loc}}(\\phi^{(n)}) = \\frac{\\bra{\\Psi_\\mathrm{T}}H\\ket{\\phi^{(n)}}}{\\braket{\\Psi_\\mathrm{T}| \\phi^{(n)}}}\n", "\\tag{5}\n", "\\end{equation}\n", "defined as the mixed expectation value of the Hamiltonian with the trial wave function $\\ket{\\Psi_\\mathrm{T}}$.\n", "\n", "\n", "The trial wavefunction can be in general a single or a multi-Slater determinant coming from VQE for example. This might help in achieving more accurate ground state energy estimates.\n" ] }, { "attachments": {}, "cell_type": "markdown", "metadata": {}, "source": [ "\n", "The implementation of AFQMC we use here is from [ipie](https://github.com/JoonhoLee-Group/ipie) that supports both CPUs and GPUs and requires the following steps:\n", "\n", "\n", "1. Preparation of the molecular Hamiltonian by performing the Cholesky decomposition\n", "\n", "2. Preparation of the trial state from the VQE wavefunction\n", "\n", "3. Executing AFQMC\n" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Preparation of the molecular Hamiltonian\n" ] }, { "cell_type": "code", "execution_count": 7, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "# Number of electrons in simulation: (12, 12)\n" ] } ], "source": [ "# AFQMC.\n", "\n", "# Generate the input Hamiltonian for ipie from the checkpoint file from pyscf.\n", "ipie_hamiltonian = gen_ipie_input_from_pyscf_chk(hartee_fock.chkfile,\n", " mcscf=True,\n", " chol_cut=1e-5)\n", "\n", "h1e, cholesky_vectors, e0 = ipie_hamiltonian\n", "\n", "num_basis = cholesky_vectors.shape[1]\n", "num_chol = cholesky_vectors.shape[0]\n", "\n", "system = Generic(nelec=molecule.nelec)\n", "\n", "afqmc_hamiltonian = HamGeneric(\n", " np.array([h1e, h1e]),\n", " cholesky_vectors.transpose((1, 2, 0)).reshape(\n", " (num_basis * num_basis, num_chol)), e0)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Preparation of the trial wave function\n" ] }, { "cell_type": "code", "execution_count": 8, "metadata": {}, "outputs": [], "source": [ "# Build the trial wavefunction from the state vector computed via VQE.\n", "wavefunction = get_coeff_wf(final_state_vector,\n", " n_active_elec=num_active_electrons,\n", " spin=spin)\n", "\n", "trial = ParticleHole(wavefunction,\n", " molecule.nelec,\n", " num_basis,\n", " num_dets_for_props=len(wavefunction[0]),\n", " verbose=False)\n", "\n", "trial.compute_trial_energy = True\n", "trial.build()\n", "trial.half_rotate(afqmc_hamiltonian)" ] }, { "attachments": {}, "cell_type": "markdown", "metadata": {}, "source": [ "### Setup of the AFQMC parameters\n" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Here we can choose the input options like the timestep $\\Delta\\tau$, the total number of walkers `num_walkers` and the total number of AFQMC iterations `num_blocks`.\n" ] }, { "cell_type": "code", "execution_count": 9, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "# random seed is 96264512\n", " Block Weight WeightFactor HybridEnergy ENumer EDenom ETotal E1Body E2Body\n", " 0 1.0000000000000000e+02 1.0000000000000000e+02 0.0000000000000000e+00 -2.2437583763935545e+04 1.0000000000000000e+02 -2.2437583763935547e+02 -3.7639365190228011e+02 1.5201781426292453e+02\n", " 1 4.2276634193515412e+02 1.4127560668989827e+03 -1.1711742028818304e+02 -2.2473358126540003e+04 9.9999999999999986e+01 -2.2473358126540006e+02 -3.7646854013277283e+02 1.5173495886737268e+02\n", " 2 1.0031922288872407e+02 3.8320523739865604e+02 -1.1743088014788954e+02 -2.2489226882493567e+04 1.0000000000000001e+02 -2.2489226882493563e+02 -3.7650504938463922e+02 1.5161278055970348e+02\n", " 3 9.9900990681040355e+01 1.0008400623205630e+02 -1.1736864885170948e+02 -2.2495677577437204e+04 9.9999999999999972e+01 -2.2495677577437212e+02 -3.7659644834889821e+02 1.5163967257452603e+02\n", " 4 1.0009188692360159e+02 1.0005173726372723e+02 -1.1748969527283802e+02 -2.2495531836556856e+04 1.0000000000000001e+02 -2.2495531836556853e+02 -3.7675907314082951e+02 1.5180375477526098e+02\n", " 5 9.9997269300807844e+01 1.0010618465796188e+02 -1.1752703012577417e+02 -2.2502732667629320e+04 1.0000000000000001e+02 -2.2502732667629317e+02 -3.7663343013337044e+02 1.5160610345707727e+02\n", " 6 1.0012131352337956e+02 1.0019003056579172e+02 -1.1770170647504112e+02 -2.2513369839216481e+04 1.0000000000000000e+02 -2.2513369839216480e+02 -3.7660812717909516e+02 1.5147442878693036e+02\n", " 7 9.9936984461419740e+01 9.9929966800671224e+01 -1.1765353928750643e+02 -2.2516138533920657e+04 1.0000000000000000e+02 -2.2516138533920659e+02 -3.7660292355465600e+02 1.5144153821544941e+02\n", " 8 9.9902337463172714e+01 9.9910800755312891e+01 -1.1761532255317621e+02 -2.2518524275281430e+04 9.9999999999999986e+01 -2.2518524275281433e+02 -3.7674246483479845e+02 1.5155722208198404e+02\n", " 9 1.0012943675389775e+02 1.0013880643723378e+02 -1.1780913595074867e+02 -2.2512465963277762e+04 1.0000000000000000e+02 -2.2512465963277762e+02 -3.7677999264623367e+02 1.5165533301345607e+02\n", " 10 9.9628730363609819e+01 9.9223106824565718e+01 -1.1749814144939067e+02 -2.2517668156221851e+04 1.0000000000000000e+02 -2.2517668156221850e+02 -3.7688306341863290e+02 1.5170638185641434e+02\n" ] } ], "source": [ "# Setup the AFQMC parameters.\n", "afqmc_msd = AFQMC.build(molecule.nelec,\n", " afqmc_hamiltonian,\n", " trial,\n", " num_walkers=100,\n", " num_steps_per_block=25,\n", " num_blocks=10,\n", " timestep=0.005,\n", " stabilize_freq=5,\n", " seed=96264512,\n", " pop_control_freq=5,\n", " verbose=False)\n", "\n", "# Run the AFQMC.\n", "afqmc_msd.run(estimator_filename='afqmc_src/estimates.0.h5')\n", "afqmc_msd.finalise(verbose=False)\n", "\n", "# Extract the energies.\n", "qmc_data = extract_observable(afqmc_msd.estimators.filename, \"energy\")" ] }, { "cell_type": "code", "execution_count": 10, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "" ] }, "execution_count": 10, "metadata": {}, "output_type": "execute_result" }, { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAlQAAAGwCAYAAABvpfsgAAAAOXRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjkuMiwgaHR0cHM6Ly9tYXRwbG90bGliLm9yZy8hTgPZAAAACXBIWXMAAA9hAAAPYQGoP6dpAACD50lEQVR4nO3deXgT5doG8HuSNOm+QVdaoCyyIwiCRRFQPKAU4chRRBFZBP2OR0ERFTeOcrSiCCoqIgLFI8oBF1RcEIEiCMgiyNpigVLoCoXua5L5/khmsrdJmzYp3L/r6tVkZjJ5M6B5eN5nnlcQRVEEERERETWYwtMDICIiImrpGFARERERNRIDKiIiIqJGYkBFRERE1EgMqIiIiIgaiQEVERERUSMxoCIiIiJqJJWnB3A10Ov1yMnJQVBQEARB8PRwiIiIyAmiKKK0tBSxsbFQKOrOQTGgagY5OTmIj4/39DCIiIioAc6dO4e4uLg6j2FA1QyCgoIAGP5AgoODPTwaIiIickZJSQni4+Pl7/G6MKBqBtI0X3BwMAMqIiKiFsaZch0WpRMRERE1EgMqIiIiokZiQEVERETUSKyhIiIi8hCdTofa2lpPD+Oqplar622J4AwGVERERM1MFEXk5eWhqKjI00O56ikUCiQkJECtVjfqPAyoiIiImpkUTEVGRsLf359Nnz1Earydm5uLtm3bNurPgQEVERFRM9LpdHIw1apVK08P56oXERGBnJwcaLVa+Pj4NPg8LEonIiJqRlLNlL+/v4dHQgDkqT6dTteo8zCgIiIi8gBO83kHd/05MKAiIiIiaiQGVERERESNxICKiIiIqJEYUJFXqqjRQhRFTw+DiIiMRo8ejZEjR9rdt2PHDgiCgMOHDwMAVq9ejeuvvx7+/v4ICgrCkCFDsHHjRovXpKamQhAEuz95eXlN/nncjQEVeZ2MgjL0fWUznvv6qKeHQkRERtOmTcPmzZtx/vx5m32rVq1C//790bt3bzz11FN4+OGHMX78eBw+fBh79+7FTTfdhDFjxuC9996zeW16ejpyc3MtfiIjI5vjI7kV+1CR1zmeW4JqrR4bDmZj3uju8PVRenpIRERNShRFVNY27rb9hvLzUTp1p1tSUhIiIiKQkpKCF154Qd5eVlaG9evX480338SePXvw1ltv4d1338Vjjz0mH/Pqq6+iqqoKTz75JMaMGYP4+Hh5X2RkJEJDQ936mTyBARV5ncoareF3rQ57ThdiaJeW9y8VIiJXVNbq0P2lTR557+OvjIC/uv5wQKVSYdKkSUhJScHzzz8vB2Hr16+HTqfDhAkT8NJLLyEwMBAPP/ywzetnz56NRYsW4csvv8SsWbPc/TE8jlN+5HXKq03/StuaVuDBkRARkbmpU6fi1KlT2L59u7xt1apVGDduHEJCQnDy5El07NjR7rp4sbGxCA4OxsmTJy22x8XFITAwUP7p0aNHk3+OpsAMFXkd87T3lhMFePlOkQ3wiOiK5uejxPFXRnjsvZ3VtWtXDBo0CCtXrsTQoUORkZGBHTt24JVXXpGPqe+GIutga8eOHQgKCpKfN2b5F09iQEVep7xaKz/OLqpEen4pukYHe3BERERNSxAEp6bdvMG0adPw2GOP4f3338eqVavQsWNHDBkyBADQuXNn7Ny5EzU1NTaBU05ODkpKSnDNNddYbE9ISLgiaqg45Udep6LGsjBzywlO+xEReYt77rkHCoUCn332GT755BNMnTpVnkWYMGECysrKsGzZMpvXLVy4EL6+vhg/fnxzD7lZtIxwmK4qFcai9DahfsguqsTWtAI8OqyTh0dFREQAEBgYiPHjx2Pu3LkoKSnB5MmT5X2JiYmYOXMm5syZg5qaGowdOxa1tbX49NNP8e677yIlJQWtWrWyOF9BQQGqqqostrVq1arFTf0xoCKvI2WoknrHYNmvp/FH1mVcKq9BeIBtkSMRETW/adOmYcWKFbjjjjsQGxtrse/tt99G79698cEHH+CFF15AVVUV1Go1tm7diptvvtnmXF26dLHZtnv3btxwww1NNv6mwCk/8jpSQNUxIhDdYoIhikBqOqf9iIi8RWJiIkRRxPfff293/9SpU7F//35UVlbizJkziI6OxgcffACdzlTSMXToUIiiaPenpQVTAAMq8kLSlJ+/Rolbuxp6UG1h+wQiohapffv2SE1NRdeuXXHo0CFPD6fJMKAiryNlqPzVStzSzRBQ/Zp+AbU6vSeHRUREDZSQkIB///vf6Nevn6eH0mQYUJHXMQVUKlwbF4pWAWqUVmuxL/OSh0dGRERkHwMq8joVxj5U/mollApBXnpmK9snEBGRl2JARV6notaUoQKAW43TfqknL3hsTERERHVhQEVep6LaVEMFAF2jDUsS5JdUOXxNS6TTi9h+8gKKK2s9PRQiImokBlTkVWp1etQYi88DjBmqQI3hd1m1tt41opxVWlWLz/dm4VJ5jVvO1xCbj+fhwZV78fqPJzw2BiIico8WEVBlZmZi2rRpSEhIgJ+fHzp27Ih58+ahpsb0ZZiamooxY8YgJiYGAQEB6NOnD9asWePwnGvXroUgCBg7dmy975+amorrrrsOGo0GnTp1QkpKihs+FdljvuyMnzFDFehrCKhE0XZZmob6375zmPvVESzbfsot52uIc5cqAQDnL1d6bAxEROQeLSKgSktLg16vx7Jly3Ds2DEsXrwYH374IZ577jn5mF27dqF379748ssvcfjwYUyZMgWTJk3Cxo0bbc6XmZmJp556CoMHD673vc+cOYNRo0Zh2LBhOHToEGbNmoWHHnoImzZtcutnJINKY8DkoxSgVhn+evr5KKEwLBNlsXByY+QVG6YPPTmNWGb8LGVu+kxEROQ5LWLpmZEjR2LkyJHy8w4dOiA9PR1Lly7FwoULAcAiuAKAmTNn4ueff8ZXX32FpKQkebtOp8P999+Pl19+GTt27EBRUVGd7/3hhx8iISEBb731FgCgW7du2LlzJxYvXowRI0a46ROSpNzY1NPPRylvEwQBARoVSqu0KK3WItIN7+MNwYwUHLorSCQiIs9pERkqe4qLixEeHu7yMa+88goiIyMxbdo0p95n9+7dGD58uMW2ESNGYPfu3Q5fU11djZKSEosfco6UoQrQWMb6Qcbn7go+So3nKa3yYEBVIwVU7pnGJCJqLrt374ZSqcSoUaMstmdmZkIQBJufiRMnWhy3evVqXH/99fD390dQUBCGDBliM6OUmpoKQRAQFhZms3jyvn375HObE0URH330EQYOHIjAwECEhoaif//+ePvtt1FRUeHGK2CrRQZUGRkZWLJkCR5++GGHx6xbtw779u3DlClT5G07d+7EihUrsHz5cqffKy8vD1FRURbboqKiUFJSgspK+7UvycnJCAkJkX/i4+Odfr+rnRQwSfVTEinAKnNTACSdx5MZqjJjIMUpPyJqaVasWIHHHnsMv/76K3Jycmz2//LLL8jNzZV/3n//fXnfU089hYcffhjjx4/H4cOHsXfvXtx0000YM2YM3nvvPZtzBQUF4euvv7Z5/7Zt29oc+8ADD2DWrFkYM2YMtm3bhkOHDuHFF1/EN998g59//tkNn9wxj075Pfvss1iwYEGdx5w4cQJdu3aVn2dnZ2PkyJG4++67MX36dLuv2bZtG6ZMmYLly5ejR48eAIDS0lI88MADWL58OVq3bu2+D2HH3Llz8eSTT8rPS0pKGFQ5SepBJd3hJ5EK090VfHjblJ8oijb/0iKiq4goArVNm0FxyMcfcOH/P2VlZfjf//6H/fv3Iy8vDykpKTZlN61atUJ0dLTNa/fs2YO33noL7777Lh577DF5+6uvvoqqqio8+eSTGDNmjMV35oMPPoiVK1diwoQJAIDKykqsXbsWjz/+OObPny8ft27dOqxZswYbNmzAmDFj5O3t27fHnXfe2eSzRR4NqGbPno3JkyfXeUyHDh3kxzk5ORg2bBgGDRqEjz76yO7x27dvx+jRo7F48WJMmjRJ3n7q1ClkZmZi9OjR8ja93nB7vkqlQnp6Ojp27GhzvujoaOTn51tsy8/PR3BwMPz8/OyOQaPRQKPR1Pm5yD6pB5V1hsq8dYI7yBkqD075SZ9FqxdRrdXD10dZzys8o6SqFqIIhPj5eHooRFeu2grgtVjPvPdzOYA6wOnD161bh65du6JLly6YOHEiZs2ahblz5zr1j8LPP/8cgYGBdmeYZs+ejUWLFuHLL7/ErFmz5O0PPPAA3nzzTWRlZaFt27b48ssv0b59e1x33XUWr1+zZg26dOliEUxJBEFASEiI05+xITwaUEVERCAiIsKpY7OzszFs2DD069cPq1atgkJhO1uZmpqKpKQkLFiwADNmzLDY17VrVxw5csRi2wsvvIDS0lK88847DjNIiYmJ+OGHHyy2bd68GYmJiU6Nm1xTYawrCnAQULmrhkoKZkq9IEMlPfbGgEqnFzFy8a/QiSJ2PnMLfJQtskqAiNxoxYoVck3UyJEjUVxcjO3bt2Po0KHyMYMGDbL4nt6xYwf69u2LkydPomPHjlCr1TbnjY2NRXBwME6ePGmxPTIyErfffjtSUlLw0ksvYeXKlZg6darN6//66y906dLFTZ/SdS3iLr/s7GwMHToU7dq1w8KFC3HhgmkJEimluG3bNiQlJWHmzJkYN24c8vLyAABqtRrh4eHw9fVFz549Lc4bGhoKABbb586di+zsbHzyyScAgEceeQTvvfcenn76aUydOhVbt27FunXr8P333zflR75qmS+MbE6qoXJXAFRaZehOXqPVo1qrg0bV/MGMebatrFqLVoHel9UsqaxFjrHFRFFFLSKCvG+MRFcEH39DpshT7+2k9PR07N27V65pUqlUGD9+PFasWGERUP3vf/9Dt27d5OfmSYv6GjTbC7amTp2KmTNnYuLEidi9ezfWr1+PHTt2WBzjrsbPDdUiAqrNmzcjIyMDGRkZiIuLs9gnXcDVq1ejoqICycnJSE5OlvcPGTIEqampTr9Xbm4usrKy5OcJCQn4/vvv8cQTT+Cdd95BXFwcPv74Y7ZMaCKmgKrpMlSiKFoEM+XVtgGVKIpYf+A8ukYHoXdcqN3z7DldiAul1Rh9bcPS9OVWAZU3sg76GFARNRFBcGnazVNWrFgBrVaL2FjT//dEUYRGo7EoKI+Pj0enTp1sXt+5c2fs3LkTNTU1NoFTTk4OSkpKcM0119i87vbbb8eMGTMwbdo0jB49Gq1atbI55pprrkFaWlpjPl6jtIj8/eTJkyGKot0fSUpKit39dQVTKSkp2LBhg80269cMHToUBw8eRHV1NU6dOlVv3Rc1nDTl5yigckfNU2WtDnqzf8jYO+eJ3FI8/cVhzFl/2OF5/vXZH3js84MNbg5qHdR5I/O2Ep6sNyMiz9Nqtfjkk0/w1ltv4dChQ/LPn3/+idjYWHz++ef1nmPChAkoKyvDsmXLbPYtXLgQvr6+GD9+vM0+lUqFSZMmITU11e50HwDcd999OHnyJL755hubfaIoori42IlP2XAtIkNFVw85Q6Wxf5efO6b8rAOD0mrbxYmlICm/1H6wVKvT42KZYemjgpJqRAX7ujQGrU6Pqlq9/Nxbm3uaB332rhMRXT02btyIy5cvY9q0aTYF3uPGjcOKFSssmnDbk5iYiJkzZ2LOnDmoqanB2LFjUVtbi08//RTvvvsuUlJS7GafAGD+/PmYM2eOw/333HMPvv76a0yYMAEvvPAC/va3vyEiIgJHjhzB4sWL8dhjjzm13FxDMaAiryJnqHyabsrPOiizl3kpMdZYlVTW2m1pYJ65kY51RbnVmoTeO+Vn+mzMUBFd3VasWIHhw4fbvVtu3LhxeOONN5xqTfD222+jd+/e+OCDD/DCCy+gqqoKarUaW7duxc033+zwdWq1us62R4Ig4LPPPsNHH32ElStX4tVXX4VKpULnzp0xadKkJi/VYUBFXsVhhsqNbROsAwN75yypNAQSetEQ/ARajUfab/3YWdaBobdmqCym/Lx0jETUPL777juH+wYMGCCX4ThTHD516lR56i4zMxNDhgzBBx98gBtvvBFKpeEf1EOHDq3zXGPHjrXZr1Ao8Mgjj+CRRx6pdwzu1iJqqOjqIdUSOayhckOtkXVgYDegMs9A2QmYzLNSDcpQOTEGb2BdlE5E5G7t27dHamoqunbtikOHDnl6OA3GDBV5lcpa+0XppqVnGl/HY71+n731/CwyUFW1iIWf1X6t3cfOsg5OvLUo3Tyb58l1D4noypaQkIB///vfnh5GozBDRV7FlKGyWhzZV2WxvzGcy1CZT+nVs79BGSrrGirvLPhmhoqIyDkMqMirVNZIa/k5yFC5pYaq1uq57TmLzTJUxXam/OrbX+8YrAIod0xlNgVn2ib8dDQXD/93f4Oug7todXqculDm8cZ+RK7g31fv4K4/BwZU5FXKjXf5OVrLr7xGC72+cX/5ncpQVdZTQ9XIonTrAKolFKWXOsjEffTraWw6lo9fT16wu785vLEpHbe+tR1b0wo8NgYiZ/n4GNbFrKjw0GLIZKGmxtACRyqGbyjWUJFXkTNUDu7yE0Wgotb2rjtXSG0TBMFwPrs1VPVM6Vnudz0Yail3+Vm0TXAwxiJjQOnJDFV6Xqnhd34pbu0W5bFxEDlDqVQiNDQUBQWGfwD4+/s7tbAwuZ9er8eFCxfg7+8PlapxIREDKvIqcobKqg+Vr48CSoUAnV5EebW2UQGVNHXVOlCDC6XVduuXLDNQrmew6h1Dtan4vqJG57X1SRaNPR0EjtLn92TRuhTMNeQGASJPkNahlYIq8hyFQoG2bds2OqhlQEVeQ6cX5e7h1hkqQRAQqFGhuLIWpVVaRAU3/H2kICEmxNcYUNXTNqHeDFXD2yZEBfvizMVyOZC0tv3kBYT7q9ErzraRXnMoq6cPlSiKchDTkOvgLnIjVg+OgcgVgiAgJiYGkZGRqK3l31tPUqvVUCgaXwHFgIq8RmWtqa7Ium0CADmgamw2RwoSooN9cRjFNsXWhiCh7hqp+jJY9TEFVBpDQGWnKP1CaTWmrNqLVoEa7Ht+uMvv4Q6l9dzlV63Vo0ZnCIIbkqlzF+nPoDHTjnq9CEEAp16oWSmVykbX7pB3YFE6eQ1p2RmFAGhUtn813bX8TKlZhsr8uaSyVgetWeF709zlZwigoo1rANoLVvKKq6AXDYFVtdYzdwGW1XOXX3Fl42rJ3ME8AG5oUFet1WH44u2YvGqfO4dGRFcRZqjIa1SY9aCylyUI0Bj+FdfYWh05QxXiZ/FcYh0g2Z/yM72mslaHGq0eajtBoCPmU37mz80VVdZYjCkyqPn/FWu5OHI9DVA9lKFyR5bs3KVKnL5QjjMXy6HXi1AomKUiItcwQ0VeQ17Hz850HwAE+hpuNW5shqrMOkNlFVBZT+HZL0q3/OJ21FLAEalmKtIYUFXU6KCzagdRVGGWBato/mBFpxflPxMAqNHqbTJlja0lc4cSN2TJpCBaFO0HjkRE9WFARV5DmvJzGFAZM1SNrqGyCqgqa3XQGjMcgG1gUF9RuuG5a2OSgrioYI28zbowvajClKG67IGAyt51tq71KvaCDFVjp18NrzNda0/WghFRy8WAiryGKUNlfyY60E3d0qUpvpgQ0/p85oGC9IUqvZ/1F2y1ViffjejomPpIWbZwfzVUxukl68ybeYbKPLhqLtJ1VqsUchsL6+lRi/YRHqqhslwmqLZBXY/dEZQR0dWNARV5jfoyVO5YfqZaq5PrbcICfOTi99Jq26mruDA/4z7L7uxSdkkQgDahfhavcZYUPAVoVPLnsgmozL7YizzwJS8FT8G+KnktxdJqx9k7T2V2zIM6rdU0pbPMg1e2XiCihmBARV5DzlA5aNoZ5Ia7/MwzLAFqU6BgHqRJX9BSQGVdV1NslsEK8fexeI3T4zALqEyZN8tA4LJZVsoTGSqpLixQo0KgdJ2sC/jNApFqrR5Vtc1/N6IzNxG4cg42ByWihmBARV6jXAqofOrJUDViakkKZAI1KigUgimYqbINmCKCNHIGy97dbMG+Pgg2Fsq7Mk0kiqL8WQM1KvnuRetAsdhiyq/5syZSEBnoq5KDWevsoHXw4olu6dZjaMiUnTfUghFRy8aAirxGpTTlp3F0l1/jp/ykL3wpkAqUp7Jsl5IJ9vVBsJ8xA2Vn7b5gPx8E+6ls9tenWquX7+gL9FU5rA3zlik/iwyVdUBlfUekB6bLrAOghmSYijnlR0SNxICKvEZ5dT1tE9xQQyUHVMYAwV6GSvpCDfbzQYif7ZSe9AUe4qcy2+/8l7D5+P19lA5rqDw95WfK5vmYiu/r69nlgcDPHWNghoqIGouNPclrSEvPBDThXX7mU36G3z4255SCp2A/HwT72mag5IDLbMrPlayGXJCuVlpOO3rZlJ8UZAb5qqAwNlq1ucuvke0j3ME6I9WQKT/zDKCn7lYkopaNARV5DSnQ8GvCDFWZ8S41qRg9yE6xtSlgUpmm/OwULRum/FwvSjcvSDf/bf65RFG0nPLzZA2VRgWlsbVDWR13+QGeye440zesPmybQESNxSk/8hqVNXVnqNxSlG5dQ6WxU0NlNuVnykDZC7jsZ7AaOgbzKb/Saq1F53SPTPmZTY/amxoFTMFH60A1AM8UpUtjCDAG4ixKJyJPYEBFXkPqFO4oQyVlkxrTNqHUesrPToaq2KIoXWWxzWK/nymD5cqXuPQ5TRkq6S4/U8sB66VmPFKUXm3bNsF6WRYpM9cmzN/w3BNF6cb3jA/3txiTs0RRZFE6ETUaAyryGlIfqgAHd/nJxds1OotGm64oc1SUXm07pRfipzJlqOppm+BaUbrl57Q35SdN8UkdyitqdDbr6DU1aTxBDjJUer0o96qKN/bs8siUn9w3zBBQuZqhqqo1La5sfj4iIlcwoCKvIQVUfj51F6UDtuveOUsOEjRWNVTG7eZBgsVdfnbaJoRY7Hd+POU2hfG2mbci49py8eF+MJYvNfsCyeYtJuw1QC2r0UKKa9uENaxjvDtIAVRcA8dQVGk5ncoaKiJqCAZU5DXqy1BpVAp53buGFqY7ylBJwUO5WZBg0YfKTtsEiz5Urkz5OQiozD+TtBhyeIBaDtqae9rP/I5I+x3lDeNRqxSICNQYtzVvdsciSxbesAyVOzqtExExoCKvUd9afoIgyIFQQ+uoSs16Kxl+WwYKUqZJrVLA10dpty2CvbsAXVl2xdFdfuafqdhYhB7qp0aov6Hgu7nv9LMsSvex2AaYT43ab4AqKa2qxX/3nMWF0mq3j9E8AG7otGOxWfAKGAL7WrMpQCIiZzCgIq8hr+Xn4C4/wHQHYEPvJrPJUFkVpZvXRwGwm4Eyb5sQqFbB2KLJ6TE5nvIzBWRS8BTqb5pWbO47/UzToz5274Y0Fe/brzWTrN17Di9uOIr3t2W4fYzFZlmyyGBfwxhd/LtRZDVl2JBzEBExoCKvUVFdd4YKML/Tr2EF2jY1VFaNPUvM7uADYBMoiKJoMeWnUAjyuZydajIVpTvuQ3VZDqjUCPOXAioP1VD5quru12WxBI9tIHL2Urnhd2G528coB7dmLSwaOuUXHqCWA0fWURGRqxhQkVcQRREVtfVnqOzdlecKuS7IQYaq2CZDZVl0Xq013REmfYHXNd1lT7nVlF+gsWaszE5Reqi/j2nKr9I2Q7UtrQBvbkpr8F2Pjuj1ot27/CprddAaP790rUL86r7bUZrqyy9x/5SfdM3NlwEqq9bKY3SGNOUXYt4ZnwEVEbmIARV5hapaPURjTFBXhipA07gpP5vFkaUArUYLvV60uIMPMAVN0pe09EWrEEzTj66u52cq9rZsm1BerYVovAjSl3yon/mUn+355317DO9vO4X9Zy879d7OMr+LMlCjksdoGKch8DWfHrV3N6SkQA6oqtw6RsC8J5ipjgtw7e+HdI7QemrBiIjqwoCKvEKF2Re41HvJnsYWpZs3qwRMU4iiCFTU6iym88x/G16rtZjmUhjvOLTXTb3uMdgvStfqRVRrDZkVaWFkw5Sf2rjN8ktepxeRXVQJADh/ucKp93aWNEYfpQCNSgG1SgGNyvC/i1LjNZQ+b7BZv66qWr1NvywpQ1VYXuP2XlrmQZ2PUiEH464ERFLmz6K4nr2oiMhFDKjIK5h6UCnlQMWeQHXD1/Or1elRVWsIWKRAyqIVQ5XW4g4+AJZf0pVaFJvV7EhcbZ1gPeVnvtSOtE8qlDZM+Und2C2n/ApKq+TlaXKMgZU9rkx/ScyXxxGMVffSNSu1KuAP8fORA13z/YBhKrfA7O6+AjdP+5mCOimj6Hrn+mKLxbBdfz0REcCAirxEfT2oJHLNUwOK0s2zWlIwY96Koay61uIOPol56wRThkpld78r45CyZEqFIGflpOm04grbgMp6yi+nyDSFllNsfzrtqz/Oo8e8TdhyIt+psUlKrWrNzMdrU8Dv6wOlWXG+9R2RNVpTQFdQ6t5pP1Mdl1TPJo2hAVN+/mqz4noGVETkGgZU5BXqW8dPEtCIonQpc+Lro4CP0vRX37y5p/nCxxLz9fys2yoY9ruW1ZCCQfPO74FmtVqiKMoZqjB/Ux8q6ym/3GJTVirXQYZqa1oBqrV6bD95wamxyWOUM1Smz2nTYqLK/vSo+dTnhTLLACqv2M0ZKqs/jxAX/ywAU8+v+orriYjqwoCKvEKllKGq4w4/wNTuoCFtE8qqbYMEw3NTMGO+8LHE/EvWbkDl61rdjXWGyvxxeY0WpdVaeSovxM8HoVKQYNWHynyazzxbZe7cZcMx5y65VmNl3V7CfIxS9sr8Lj/ANCVoHoxYT/Hlubkw3SaoczFbCJhnqOouriciqgsDKvIKUpDhbIaqIXf5mbcBMGfeY6muDJRhys9UiG3a7/w0kU4vorLWsg+V4bFSHoM03efrY+jWLk/5VdY15Wc/Q3XeGEhJgZWzrBugArDpll5iVU8mXadSiwyVZUDl7jv9zLu1m/92JUNVZBYYmrKNLEonItcwoCKvIAcZ9WSoGnOXX1mVbWbI/HlptdambYL545JKrUUhtu3++r/EzdsRmNeLBZgV28td0v3UFr8ranQWd8mZT/mVVmnlNe0kFTVaFJYbslrnL1fILRmcUWonixbsazndauoB5Tg7ZJ2hcn9A5aBvmJMBlV4vWvyZsg8VETUUAyryCtIUXn0ZqiCrwmhX2AsSACDQ15R5sW6bAJgCCYuidHtTfk5kzcot2hGYPmugWS8qU8sE01SadONjsVkdlfU0X65VYfp5s6xUVa0eF8ucX7rGbobKURNUm4Jw0xilDFVkkGHx5DwHxfMNZX2TgKt9pMrM1gKsb01CIqK6MKAiryD1oQpwcsqvURkqX/sZKos+U77mU3rmNVR27gKUprqcyGpIYwiwCurMl58xb5kAAAqFYGruafYeUoZKbSywt26dYF03dc6FXlVSFspRDVWtTi/fmSlnh+xmqAwBVK82IQCaLkNl3YjV2Sk7KUDVGBfDdrVJKxGRpEUEVJmZmZg2bRoSEhLg5+eHjh07Yt68eaipMf2LOzU1FWPGjEFMTAwCAgLQp08frFmzxuE5165dC0EQMHbs2Drf+6uvvsJtt92GiIgIBAcHIzExEZs2bXLXRyMjuQ9VfVN+dhbpdZa9IAEw1VAVV9bKmS97bROKK2vtF62b3QVY/xikwNF+QFVerZOLz6WpPgCmO/2MU3hVtTo549SzTTAA24zVeau6KevndbHuKA9YZqjM66Sk62eaLrOtoeoVZwio8kqqXJp6rI/1UkGuBkTFVsGr6c+aNVRE5JoWEVClpaVBr9dj2bJlOHbsGBYvXowPP/wQzz33nHzMrl270Lt3b3z55Zc4fPgwpkyZgkmTJmHjxo0258vMzMRTTz2FwYMH1/vev/76K2677Tb88MMPOHDgAIYNG4bRo0fj4MGDbv2MVzu5D1U9GSo5m9SQovR6MlS5xZXy8jdBvrYBU4mjtgpmmZn6goVyOy0TDM+NfahqtHJ7hLAA03tYF6ZLU2e+Pgp0iwmWx2/OJkPlwp1+9vpQBdm5GzJQo4LKmCGzN10m1VBJGaqqWr3THeXro9XpUS5lyaxaNzhblF5UYZXhYh8qImqgutMBXmLkyJEYOXKk/LxDhw5IT0/H0qVLsXDhQgCwCK4AYObMmfj555/x1VdfISkpSd6u0+lw//334+WXX8aOHTtQVFRU53u//fbbFs9fe+01fPPNN/juu+/Qt29fu6+prq5GdbWpGLekpMSZj3lVk6b86lrHDzB9wVfW6qDTi1A66Kp+NLsYv2VcxEODO8jHOKyhMj7PLjIFKeb1TXbbJtiZ8qvViaiq1ddZB2a9OLNpDOYL+0otE8wyVHLrBMP7S3f1xYb6ITbUz7DNKkMlTfGF+fvgckWtSxkqewX85r2yTMXg9ttLSKQMVXy4P0L8fFBcWYv8kiqLov6GMs+SSeNwte2BdesH6c+yRqtHVa0OvnUsg0REZK5FZKjsKS4uRnh4uMvHvPLKK4iMjMS0adMa9L56vR6lpaV1vndycjJCQkLkn/j4+Aa919VEylD5a+qO8c3vjKurMH3uV0eQ/GMavj+SazreUYbK+FyqQbL+sjf/krZ3F2CAWikHbfV9kVsvO2P9ucqrtfLaclJWyvBYau5ZYxyrIXiKDfFDbKivxfglUgB1Q4dWxueu1FBJLSbMGntqTG0RrPs/GR6bMnkAUK3VyRmgyCANooMN43RXYboUDAWolaYsmYuNOU0BleH6BqpNNwAwS0VErmiRAVVGRgaWLFmChx9+2OEx69atw759+zBlyhR5286dO7FixQosX768we+9cOFClJWV4Z577nF4zNy5c1FcXCz/nDt3rsHvd7VwNkOlUSnlImxHhell1VocyykGAOzPvGSxHbAMEgBTdkNaxDfYer/ZNJK9PlWCIDh9u73UNiHQaokd87v8pCAkzN/xlJ/UGT021BexIYYMlaMpv8SOrSyeO0MKPoMs+lDZa4Bqr8GpYZ90PdVKBUL8fBAZbLjTz12F6XUGdZVap2q1zBdGBgw3AASxWzoRNYBHA6pnn30WgiDU+ZOWlmbxmuzsbIwcORJ33303pk+fbve827Ztw5QpU7B8+XL06NEDAFBaWooHHngAy5cvR+vWrRs03s8++wwvv/wy1q1bh8jISIfHaTQaBAcHW/xQ3eQMVT1F6YBZE0wHAdXh80XyrfB/ZF2Wt9vr/g3Ydk4P9rMOuAzPL5bVQGs8sXlRuvlr6stq1FeUXlqlRZG8FIr5lJ/hcZE85WcISmJCzKb8ik0F38WVpmyalKHKLqqEXu9cQXiZnelRywao9haJtrwGUkAVEaSBIAhyhsptAZWdMUiBUY3ZQth1sZ7yA8xvMmBhOhE5z6M1VLNnz8bkyZPrPKZDhw7y45ycHAwbNgyDBg3CRx99ZPf47du3Y/To0Vi8eDEmTZokbz916hQyMzMxevRoeZteb/gfrkqlQnp6Ojp27OhwHGvXrsVDDz2E9evXY/jw4c58PHJBRbUUUNVfsxLoq8LlilqHAdXBrCL58YncUlTUaOGvVtm9c006n7lg6+fGL1hpORiV2WLGptc4VwztaMrPfOkZKbgMtZOhKq6UpvxMGaqoYF8IgqHup7C8Bq0DNfL0XqsANTq0DoBSIaBWJyK/tAoxxoxWXaQmofYWRy4168cVYjdDZfiMBcaAqrWxB1V0iHHKz00Blb1gKMA4ZacXDYFdfX3NzBehlhg+RyWn/IjIJR4NqCIiIhAREeHUsdnZ2Rg2bBj69euHVatWQaGwTa6lpqYiKSkJCxYswIwZMyz2de3aFUeOHLHY9sILL6C0tBTvvPNOnXVOn3/+OaZOnYq1a9di1KhRTo2XXFNR69yUH2DWVdzB3WJ/nDVlpXR6EUfOF2Ngh1Z1FITbzzbVtV8QBKttti0D7HF0l59F2wSzhZEl0hf+5XLjlJ9ZUbpapUBEoAYFpdXILapC60ANzl0y7I8L94dKqUBsqC/OXarEuUuV9QZUoijaX8tP6lJfo5Nruey1j6is1aFGq5czVFJTzyi5hso9CyRbN/UEDFN2wX4+KKowTM9K7+mIvaCMvaiIqCFaxF1+2dnZGDp0KNq1a4eFCxfiwoUL8r7o6GgAhmm+pKQkzJw5E+PGjUNeXh4AQK1WIzw8HL6+vujZs6fFeUNDQwHAYvvcuXORnZ2NTz75BIBhmu/BBx/EO++8g4EDB8rn9fPzQ0hISJN95quNKUNV/19JaerJXg2VKIo4eK4IABAX5ofzlyvxR1aRIaBykKGyXtvPuoZKpVQgUKMy9ajytR2js4vyljps7KmU9xdVOC5KN9VQmab8ACAm1A8FpdXILqpEr7gQOUMVF2bYHx/mj3OXKnH+cgUGJNR9M0dlrU6eMrWXoTJ//2CLonXT/tKqWjlDFWEVUBWUumvKz7aeTXpeVFHrVOsE6z5U5udjQEVErmgRRembN29GRkYGtmzZgri4OMTExMg/ktWrV6OiogLJyckW+++66y6X3is3NxdZWVny848++gharRaPPvqoxXlnzpzpts9H5jVUTkz51dHcM7OwApfKa6BWKTBhQFsApjoqR4sjWwdY9m7pt5fBsLe/3qL06rqL0gvLqy2WQpGY2ibUoKSqVv7s0h1+bYy/pcyVdIdffJg/AFNgJWWu6iIFngoBFlObvj5Km67s5mOUAk/AcKefdYaqqe7ys84outI6QapJq+tuRSIiZ7SIDNXkyZPrrbVKSUlBSkqKS+e1d7z1ttTUVJfOSQ1T7uRdfkDdy89I03292oTghg6GTMzBrMvQ60W7hdbSewoC5Kae1gXngHWjT9uAylSQXc+UX42jDJXhuTQGXx+FRQ8kecqvolbODoX6+8gZvRj5Tj/DPumOvvhwU4YKcK51gnm/LuupzUBfFS6V1yDbGFDZFvCr5D5VF4yZKDlDFWL4fbGsGlqdXm510FD27vIzPHe+c73donRmqIioAVpEhoqubKIoolLqlF5PHyrA8m4za1I26rq2oegRGwIfpYCLZTU4kWdqrmpdQyUIgkWQZT2FBNhvD2Cx38m2CY6COuvn5vVTgGnKr7JWhzMXywHAohYqxljwLQU6UoYqTspQGQMrZ9bzM7VMsP2c0jilO/VsC/hNvapMGSrD2FoHaKBUCNCLpoafjWG6y89+hrG+ejbAbMrPXg0Vi9KJyAUMqMjjanR6uR1BfXdlAWZF6fYyVMY7/K5rGwZfHyV6xBrq3H49eRGAoSeSeRd0iXnxtd0MlK/tlJDFfie/hMsdBFQalQIqs67v1tOKQRpTw8kTuYbgUJrmMzw2ZqiKKiGKohw4xYdZZqicmvJzMEbzbfamJQHLWrICqyk/hUKQH+eXND6gspddMh9DfRkqrU4vf1bLtgmuLV9DRAQwoCIvIGWnAMDfiaU+zJdAMVdWrUW6MRN1Xbsww++2ht+/nrxg8VpH5wQcZajq2e902wT7mThBECy2mRdJA4ZgRPrSP24MqCwyVGbLz1wqr5Fr0qQeVfHhhoAqr6QKWl3d/ZnstUyQ2LSYcDDdVlRRi4tllkXpgPmdfo2vo3I05edsPZv59Ky9PlTOZLiIiCQMqMjjpAVu1SqFU3U15h27zR0+Z2jo2SbUT/7ivq5dKABg/9lLFq91dE7AQQbKTgNLi/1OfgmXOehDZT0G6yk/823Hc4wBlVmGKjbEdAddZqEhOxUVrJHrsCICNVCrFNDpRbnOyhFH/boA26aojrJDWZcqUGtck7B1oCmgcmdzT4d3+TmZYZLupjRf4Nn8fJzyIyJXMKAij6t0oSAdsFymxZxUP9W3bai8TcpQSV/uDgMqO922zdmbErK3v64vYVEUHU75AZbrFFpnqAAgxLhNqpOSpvkAQ9DiozTUJx0wBo/SNB9gyHDFhUp3+tVdR+WoX5e9bbYZKsPzjIIyAIblc9Qq0/9m3NncU+pkbhvUSXfp1R0QOZwyZB+qK9Khc0XIKnR++SUiVzGgIo+Tp8Gc6EEFWC7TYk6qn+prDKIAw5RXtFlzR0dTfkEuFaXX0Yeqji/haq2pVixAYxs8mmetzJedkYRaffGbT/kpFIIcrOw9YwiopFYJkrhw6U6/uuuo5KL0erJoSoWAAKsgWLph4PQFQ0AlFaRL3Lmen73GnobntkXp5dVaJP9wQl7jETD19HK0GDZrqK4c+SVVGLd0Fyau+N3TQ6ErGAMq8jip3seZgnTAvGO36QtTFEUcNLvDz5w07QfYDxIA+2vWmQt2oW2Co0V5zaco7QWPllN+tu9hPQ0YG2oZrEiLJEsBlVQ3JZF7UdVzp1+dRekWtWa2bRWkwPKsMQtmXj8FuG/Kr8rYjR1wnCUzD4iW7ziNZb+exkvfHJO3ldRT1F7XnyW1LH/ll0GnF5F1qUJeiJ3I3RhQkcdJ/4OzznY4ItdQmWWoMgsrcLmiFmqVQr6zT9I33pSxqq8oPUCttFvHVX/bBMM2nV6UA0Rr0nSfv1oJhUKw2W8eZNU15QcAggCbZVWkAnSp2No6Q2XqRVV3hqq02nHbhHrvhrRa9zDSQUDV2KJ0KRgSBCDQKji1nn4VRRHfHsoBYJgWlto52OuSbv0ZHP1ZUsti3n9NakpL5G4MqMjjXM5QyUXppi8784ae5jU7gGWGqr6idHtBAmAZRIXYKVr39VHAR2kIkhxNFdWV+QEsgz37U36mbZFBGvhYBX5SLyqJeQ0VYGryWW8NVVUdNVT1TY1abbPOUEWFSBmqxrVNkKf7fH1sglPrOy6PZpfgtLF3lygC29IKAJi6pFtnqPx8lHILi8YUph/LKcbNb2zDhoPZDT4HuYf5PyLq+wcFUUMxoCKPM2WonKuhMgVUpi+7PxxM9wGQG3wCddRQGbfbCxKA+tsmCIJQ791hjhZGltQ75Rdg2hYb6mez33qb7ZSfsReVk1N+9qZHg+op3rcOSG0CKmOGqqxaK79PdlElhi/ajue+tly8vC5SQbq9OzKlcZVVa6HXi/j2T0NAI/0d+Pl4vvEc9gMqQRDcUke16Wgesi5V4CsGVB6XbZaVymaGippIi1h6huwrKK3CgFe3eHoYbuNqhqqqVo/2z35vse86s4J0idTg89C5onprqOwFCdbbHWWxQvx8UFheg7OFFegaHWyzv7yOlgmG7eZ3+dlmqMzHEBtiL6AyZagUgumOOonU5DO/pBrVWp3dBqdAPRkqi1qyuttLALYBVaBGJS80nV9SBf9WAZiz/k9kFJQho6AM025KQMeIQLvjMueoZYL5uETREBB9+6dhuu9fwzpj8S8nsTPjAiprdKYMlZ3gNdj4Z9mYXlRSVuzMxbIGn4Pcg1N+1ByYoSKvIAjAwA6tnDo2xM8H3WNsA5bWgRokdrR/jjF9YqFWKeSGn9b6tA2Fn4/S4eujgn3RISIA17UNtVhjz5x0d+EzXx5Gel6pzX5TDyr7r6+rsadhmynIsp7eM2zzs3hsPSUYHqCWW1Nk1zHtUVrH1GT9d0Navsb6Lj/A0B8LAPKLq/DfPWex61ShvO+/u886HJc5KQtoLwDWqJTw9TF89s0n8pFfUo1gXxUeGdoBbUL9UFWrx86Miw4zVIbP5txSQnWRlgg6f7kS1VrWYnmS+TRfXX/3iRqDGaoWrHWABgdeGN7s79sU9z35KBUOs0PWFAoB3z12Ey4bGzNKgn19bOqnJFNuTMDEG9rZBBmSrtHBOPzvvznc76NU4OdZN0Mh2BaTS14e0wOnLpTh0LkiTFzxO9Y/nIj2rQPk/XX1oLLebu9amLdNqG/KT6qXMicIAuLC/HAyvwznLleig4NMkDSVWl+Gyn4gYrlNapNgLirYF6culGPP6UJ8tOM0AOCOXtH44UgevjxwHnNGdKl3Tce6MlTS9qraajlAu6NXDDQqJW7rHoWUXZn45Xi+w7v8AOeXEnJEFEU5oBJFIKuwAp2jghp0LmqcGq3eou8Zp/yoqTCgasEUCgGtAm2/sK4GSoVg0YHbGY6CJWf319fFPVCjQsqU63HvR3uQlleK+z/+HV/8X6KcOaqrSzpgqiHz81HazYKZt02wbpkAGLIqAWolymt0cr2Utfgwf5zML7OYArHmbB8qe1Of1i0nrKf8ANOdfu9ty4BeBG7s1ApLJlyHE7nbceZiOb4+mI2JN7RzOD7AVNtkb9oRMARJBaXVOJJt6Dt1Z59YAJADqi1p+XLGL9TODQCNXc8vv6Ta4g7B0xfLGVB5SF5xFcy7XzBDRU2FU35EbhTqr8Yn0wYgoXUAsosqMX7ZHsxZ/yfmrP8TGw4ZipMd11CpjOdwUKPlX3eGShAE09p9jgIqY6H6qYJyu/tFUTQtPVNfDZWd/SqlQm5/4eujsBuUSXf66UVD0PbGP66FUiHIQdR/d5+tt/+T1BrCUVbTPNiLCtZgYIJhKndAQjiCfFW4WFYjd3OvK9PW0Bqq01Z1U6cv2L/e5y5V4FJ5jd195B7SPx6kYN+Z9SyJGoIBFZGbRQb54tOHBiI2xBdZlyqw/sB5rD9wHkezDWvwRQfbZpcAU11UGzvBEmAIPoJ8VVApBLQNtx8wSQXdnaPsT+d1iTZkSVb+dgZPf/GnvJ6dxLybu/0aqvqL86XtEUEam8afgOXnf2l0d/nz/qNfHPx8lEjPL8XvxuakjtQ/5Wca++jesVAa2yD4KBUY1iXS4lh7Aay8NqPZlF92USV+P11oc6w90nSf6bltYXp+SRVGvP0rJn7M7t1NSaqf6hMfCrVSAb3onqWPiKxxyo+oCbQJ9cOGR2/Ed4dz5Y7egKEgfWzfNnZf0zsuBO/fdx26xdifGlIoBKyeOgAV1Tq7dwEChgBlRM8o/K17lN39/+gXh8Pni/H53iys238eW9MK8GJSd4zqFQNBECymuOy1sfD1UUCpEKDTi3X27MotrrJbkA4A/duHQSEY6pr+0S9O3h7i54Oxfdvg871Z+GR3Jm4w3qRQXq3FR7+exvnLlRjSJQK3dI00m/Kr/67MMX0sr/fw7lHynX+OzmG9lFBVrQ7jl+3G+cuVeGNcb9xzfbzd95WcMWakIoM0KCittgmwAGB/5mVU1OhwPLcEF0qr7U6PUuOdN9ZMxYX5Iyu0AmcLK5B9udLhtDhRQzGgImoikcG+mHZTgtPHC4KAUb1j6jzGXlsIc7Ghfvh73ziH+32UCiTf1Qt3XdcGz311BH8VlGHm2kOYufaQxXGBGpXdbu6CICBQo0JxZW0d022G/61EOKhx6xEbgkPz/oYgje3SNZMS2+HzvVnYdCwfecVVOJ5bjBc3HJMLib/84zzUKgXUxnq2+qb8OrQOQM82lneEDu0SAZVCgFYvQhDs14pZ96FK2ZUpZzpe/OYoeseH2G2NIZECqFu6RmLtvnN2A6qjZusKHs0ptsmckXtIU35xYX5oE+pnCKhYmE5NgFN+RFeh69uH4/vHB2P2bdfILQbM3dSptcPXDu7cGlHBGnSOtD+tKGV37N3hZ36MvenAbjHBGNA+HDq9iHuW7cbUlP3ILqpEm1A/TL0xAQmtA1Cj1csF/uEB9jN1fY0NXiff2N7ueoNS9ivEz7bTOmB5l9+l8hq8vy0DABAb4otqrR7/XPOHxdqM1qQAalhXQ5B0sazGpsD9aLYpoDpm9pjcSwqE48L85BpD9qKipsAMFdFVSq1S4LFbO2PGkA6otFqzrq4WFksm9IVOLzq861GajrRea9BZkwa1w97MS8i6VAGlQsC0mxIwa3hn+KtVeDGpG9LzS/HjkTzU6PS40UHgN7ZPG9zUKcLhNNpt3aOwM+Oi4wyX3IdKi3e3/IXSKi26xwRj9dQBGL1kJ05fKMfzXx/B2+P72ARstTo9sozL+/SOC5Gn/TIvluPa+FAAhuL/Yzkl8muk+jpyv2yzgEqq12OGipoCAyqiq5xGpXTYNd0eQRCgUjruxzXlxvYAgLuus18rVp8RPaJxS9dIVGt1eO6ObhaLXQuCgK7RwXVOt0nH1VWTNKZPLH44kovh3ezXmkkZqnOXKnAy39Ck9flR3RARpMGS+/ri3o/24JtDObihQytMGNDW4rXnL1dCqxfh56NEVJChIaxURyUFVLnFVRZ395lP/5H7aHWmHlRxYf5oY1wtgOv5UVNgQEVEbtWzTQjeuufaBr/eR6nAysnXu3FEtkL91fjfw4kO90vTllLX+KFdIuRs2PXtw/HU37pgwU9pmPftMQxICLdYLke6o6996wAoFAISWgdiz+lL8lI0gGm6Ly7MD+cvV+L85UoUVdQ4vNmAGia3uAo6vQi1UoGIQA3imKGiJsQaKiIiK+ZTgQoBmHt7N4v9D9/cAYM6tkKNVo8NVosfSz2nOhi75Eu/zQvTjxqn+27o0ArtWhnuNjOfAiT3kDJRbcL8oFAIFjVU9fU6I3IVAyoiIivmHd/HXx8v9++SKBQC7rrOcDfl9pMXLPZJmagEYyCVIAdUpl5UUhF6z9hg9DROaR5lYbrbSXf4SbVTMcYVBqpq9WyoSm7HgIqIyIqvjxLdYoIREaTBE8OvsXvMzZ0NU4BHsotRWFYtb5d6UMkBVUSAvF3Kikg1Uz3bhKCHsa3DUTsZqqPZxfjxSK47PtJVKbvIVJAOGOoFI421dZz2I3djQEVEZMc3j96IrbOHINLB3YqRwb7oFhMMUQR2ZlyUt0tTe1IgFR/mD6VCQHmNDgWl1SgorUJ+STUEwdAmQspQWbdOqNHq8eDKvfi/NX843aGdLJm3TJBI035c04/cjQEVEZEdapUCQQ6WtpEMuSYCALA93TDtV16tle8qk2qn1CoF4o1f6KcvlOOYsUVCx4hABGhU6BFryFCdvliOUrOlblLTC1BonJb68o/z7vpYVxV5ys8soJIeM0NF7saAioiogaSA6te/LkCvF5FZaMhOhQeoLe7YSzArTD9qVj8FAK0CNYg1ruN4IrdUfs3XZsXuPx7JQ1WtZa8wql+22bIzEt7pR02FARURUQP1axeGALUSF8tqcDy3xDTdZwygJAmtDW0Vzlwss6ifkvRoY1mYXlxRiy0nCgAYCuRLq7XYfDzf5v0LSqrw8Y7TKK+ja/vVSqvTI7dI6kFlJ0NlNeW3+Xg+lmz5C3o97/6jhmFARUTUQGqVAoOM/am2n7xgU5AukQvTL5bLXdHNG5bKd/oZg62NR3JQo9Oja3QQJg9qDwD4ymraTxRF/Ouzg/jP9yfw72+PufmTtXz5pdXQ6kX4KAWLhbpjQ2wzVGXVWsxcexBvbT6JX07YBq5EzmBARUTUCOZ1VNYtEyRSPdWhc8XyF3n3WFO3915xxjv9jBmqr/8wTPfddV0b/L2voeP8r39dxIVS092Em47lYW/mJQDAF3+c9+q2Czv/uohzxuV4GmJrWj7+ueaAS2vwnTe+X0yIH5Rm6zVKGSrzc238MwcVxuWXvv0zp8HjbKzsokr8cjzfoz2yqrU6ZBSU1n8g2WBARUTUCFJAdSDrMg6fLwJgCqAkUoB10dheoV0rf4vmoVKGKqOgDOl5pdh/9jIUAjCmTxt0iAjEtfGh0OlFfGf8sq/R6pH8YxoAINTfB6II/Of74w36Ii6tqq1zoefG+vrgeUxc8Tv+/sEu5BsL9l2RU1SJxz8/hB+O5OFfn/0BrU7v1OusWyZIpIDqckUtKmoMn3vd/nPy/i0nCuTtzUkURTy0ej8e+mS/R4O65746iuGLfsVPR9muw1UMqIiIGiE+3B8dIgKg04s4dcGyZYIkOtgXvj6m/932NJvuAwwtGCKCNNCLwGs/nAAA3NiptbzA9DjjuohSofonuzNxtrACEUEarH84ERqVAntOX8LPduqs6nLkfDFufH0rhr6ZiqxC+xmkkqraOrNfNVo9dp26aLdovqCkCv/+9jgAQzD56Jo/UGsVEFXUaDFn/Z+Y8cl+FFfUWuwTRRFzvzoiB3x/ZBXh3S1/OfXZ7LVMAAzLCkmNW7MvVyKjoBR/ZBVBqRAQHeyLylqd3Xq1prb/7GWcyDVMB3+w7ZRHslR5xVXYcMjwd2zlb5nN/v4tHQMqIqJGurlzhMXz9q0sAyppTT+J1MzTnHTXn9R53Xxx6aTesVApBBzJLsbeM5fkoOKpv12DzlFBmD64AwAg+YcTqNGaApZanR7b0gtQYCczlFFQhgdX7UVJlRYXy6ox47/7bYrbz12qwO1v70DSkp14+5eTNueoqtVhaso+3Lf8d0xasdfi9aIo4rmvj6C4shZdooIQpFFh/9nLcsAIAEUVNZj48e9Yf+A8fj6ejykpluf44sB5bD95AWqVAk/eZmiwumRbBvY40ZdLaplgfoefROqcfr6oEv/bZ8hODesSibv7G7rff+eBDNGne87Kj9PzS7EtvaDZx/DZ3izojEX5e89cwukLZfW8gswxoCIiaqQhXUwBVZtQP/j6KG2OMZ8G7NUmxGa/+V1//molRvSIlp+HB6gxtEskAGDGf/ejpEqLrtFB+Ee/eADA/w3tiIggDTILK/DJ7kzo9CI2HMzG8EXbMWXVPgxbmIpVv52RvyzPX67AAyt+x6XyGvRsY+gIn5ZXitnr/pTvcjt/uQITlu+Rp87e/uUvi+xQVa0OD//3gNzUdG/mJUxJ2SdPl31zKAe/nCiAj1LAOxP6yAtmr/otE9/+mYO84ircs2w3/sgqQrCvCsG+KvyRVYQZ/92Pqlod8kuqMH+jIbv1xPBr8PitnXF3vziIIvDE/w6hqKLupWPkdfxC/Wz2SdvOXizHV8Z6tfHXx+POa2MBGIJa6/NfLq/B+v3nmmR6tLCsGj8eyQMADDZ24F+aesrt71OXGq0en/2eBQAI8zdMR6/bz/5nrmBARUTUSDcktIJaZfjfaQer6T6JeaF6j1jbgMp828ge0fBXqyz2SxmrIuO02AujusvF1gEaFZ76myGD886Wv3D7O79i1v8O4WxhBdRKBcprdHj5u+P4+we/YfvJC5j48e/ILa5C58hAfDJ1ID6c2A9qpQI/HcvDkq0ZyCmqxITle3D+ciXat/LHo8M6AgAWbT6J97dloEarxz/X/IHtJy/Az0eJl+/sgSCNCnvPXMLUlH04W1iOecY7D2fe2hldo4Pxtx7R+OdQw3me+eIwxi3dhZP5ZYgK1mD9I4OweuoABKiV+C2jEP/67CCe++oISqq0uDYuBNMHJwAA/n1nD3RoHYDc4io88+XhOqfFHNVQAaY6qjW/Z6GwvAYRQRoM6xKBzlFB6BodhFqdiE3H8uTja7R6PLhqL+Z8cRgzPtlvM23ZWOsPnEeNTo/ecSFYePe1UCsV2Jd5GfuNNx00h5+O5eFiWTUigzSYP7YnAEOG0N2f9UrGgIqIqJH81EoMTAgHYHuHn0Ta3ibUD+EBapv9Pc2mAaWFl83d0jVSrv25pWskbjJmMiT/6BeP7jHBKK3S4mR+GYJ9VZgzogsOvDgcr/69J4J8VTh8vhgPrtyLzMIKtAn1w3+nDUR4gBr92oXhP8Yv0cW/nMTY93/DuUuVaNfKH5/PuAFzRnTF0yO7AADe3JSOUe/uwNa0Avj6KLBicn88OKg9Ppk2AIEaFfacvoQRb/+K4spa9GwTjIeHdJTHOPtvXXBTp9aorNUhu6gSCa0D8MUjg9AlOgh924bh4wevh0alwC8n8rElzZDdeuMf10KlNHxVBWhUeHdCX/goBWw6lo953x5DZY1t7ZZOL8p38cWFO57y+6ugzHi928jvMdqYpTIvDF+0+SQOnzfUke06VYh/f3vMqRqnyhpdvUX0er0oZ4YmDmyHqGBfOXj+cHvzZak+2ZUJALhvYFuM6BGN1oFqXCyrxta05p96bKkYUBERucH/De2ILlFBGGcnGAKAW7tFol+7MDnbYq1NqB/uuq4N7ugVjcSOrWz2+/oo8cTwa9A9JhgvJnW32a9UCFgwrjf6tg3Fo8M6Ysczt+DRYZ0Q5OuD+we2w5YnhyCpdwwAoHWgBmseGojoEFN/pnuuj5d7XhWUViM+3A+fT78BMca+Tf8c2glzRhiCqr8KyqBWKfDxpOsxqKMhsOvbNgyrp16PALUSVbV6+CgFvPmPa+GjVFiM8Z17++DauBDc1Kk11j+SiHizgCexYyssnXgdVMbM28xbO6NLdJDF5+zZJgTP39ENAPDJ7rMY9e4OHMy6bHFMQWkVanUilAoBUcbFkM3FWk0D3tM/Xn4sTfvtPlWIgtIq/JZxEct+NQQ2kwe1hyAYMluf7D4LRwpKq/Dvb4/h2ld+xuj3frNbwybZkXERWZcqEOSrQtK1hj+fGTd3gCAAv5woQHpe07cwOJZTjP1nL0OlEHDfgLbwUSowrp/h77FUY0b1E0RPNry4SpSUlCAkJATFxcUIDrYtRiUiai7Hc0oQFaxBq0DbQEOr0+PpLw4js7Ac707oa7ege/mvp7H+wDk8P6q73DLC3P7MS/jP9ydw34C2uOf6eJv9ztiXeQkncktw34C2cubI2rb0AjzzxWEUlFZDIQAzbu6IuDA/HMspwaFzRTiRW4K4MD/sfOYWm9f+kXUZd32wCwBwffswrH9kkMX+v3/wGw5mFeHxWzvjf/uykF9SjQkD2iL5rl74cPspvP5jGpQKASlTrsdgsxsSCsuq8eH2U/jvnrOoqjVlptqG+2PNQwMtgkfJ9E/2Y/PxfEwe1B7/vrOHvP2faw7ghyN5uOu6NlgwrjcOnL2MbemGAOvWrpEYf31beZpZIooiMgsrEOrngzA7WVBHnv3yMNbuO4ek3jF4777rAACnL5Thlre242nVWjwSkArFoMeAIU87fc4rhSvf3wyomgEDKiIi9yuqqMG/vz2GDYfs35U38Ya2+M/YXjbbC0qqMOC1LQCAN//RG3f3twz8Vv12Bi9/d1x+3jEiAN89dhP81SqIoojZ6//EV39kI8hXhVu6RiKnqBI5RVXIK6mSC//7tg3F5EHt8dbPJ5F1qQJRwRp8Om0gOkeZMm65xZW48fWt0IvA5idutth3+HwR7nzvNygVAvzVSpRWWRbDtw33x5O3XYM7r41FWY0W3xzMxmd7z+FEbgnUSgVG9IzG/QPbYmBCOARBwOkLZfjxaB5+OpqHWp0ed/SKwd/7tkGwrw8GJv+Cqlo91j+SiOvbh8vvcc+y3bj13BI8rPoeSPwXMOJVZ/9orhiufH+r6txLRETkpUL91Xj73r4Y0SMay3ecRrCfD7rHBKN7bDC6xwQ7rGdrHahB1+ggVNXqcEevGJv9o3rHYP7G49CLgFqpwDv39pVvEhAEAcl39ULmxXL8kVWEb6yCud5xIXjitmsw9JoICIKAxA6tMHHF7ziZX4a7l+3Gq2N7IchXBZ1exE9H86AXgYEJ4RbBlOE8obipU2vszLiI0iotwgPUGHJNBBJaB+CT3WeRdakCs/53CO9s+Qt5xVWoNPYBUyoE1Oj0+O7PHHz3Zw46RgTAR6lAmtXUYVpeKRZtPon4cD9U1RqWOerfLszimHuvj8fZLMO0sFhTAQFUF2aomgEzVERE3kWnF1Gr09ttcQEAU1btxbb0C3hhVDc8ZOzzZe5SeQ1SdmUiUKNEbKgfYkP90CbUD5FBGgiCZehxubwGk1P24c9zRXbfa8mEvnIxvLmCkir8dCwPveNC0atNiHxXZ0WNFqt+y8SH20/JmavOkYGYMKAt7rquDc5frsSa37PwzaFseUkdlULAoE6tcUfPaKiUCnx98Dx2nSqEFAEk39ULEwa0tXj/yhod3n/tcTyFT1GQMBaRD652fEGvUFfclF9mZibmz5+PrVu3Ii8vD7GxsZg4cSKef/55qNWGeeLU1FQsXrwYe/fuRUlJCTp37ow5c+bg/vvvt3vOtWvXYsKECRgzZgw2bNjg1Dh+++03DBkyBD179sShQ4ecHj8DKiKilqWooganLpThurZhNgFSQ5RVazHvm2M4kl0EpUIBlUKAUiGgU2Qgku/qZVG878oYNx/PR/vWAejfznacpVW12HQsHwrBcGdoqL9lXVVucSU2HMxBRY0Wj93S2aYmCwC+/Xg+7jy/EIcCB6PPUxtdHmNLd8VN+aWlpUGv12PZsmXo1KkTjh49iunTp6O8vBwLFy4EAOzatQu9e/fGM888g6ioKGzcuBGTJk1CSEgIkpKSLM6XmZmJp556CoMHD3Z6DEVFRZg0aRJuvfVW5OdzNXIioitZqL8a/dqF13+gkwI1Krm5qbuE+qtt6r/MBfn64B/97N91ChgWjv6/oR0d7geA/p3bAOeBruH2M3lk0iIyVPa8+eabWLp0KU6fPu3wmFGjRiEqKgorV66Ut+l0Otx8882YOnUqduzYgaKiIqcyVPfeey86d+4MpVKJDRs2MENFRERXvuPfAOsmAfE3ANM2eXo0zc6V7+8W24equLgY4eF1/+vB3jGvvPIKIiMjMW3aNKffa9WqVTh9+jTmzZvn1PHV1dUoKSmx+CEiImpxfIyF/bXlnh1HC9AipvysZWRkYMmSJfJ0nz3r1q3Dvn37sGzZMnnbzp07sWLFCpeyS3/99ReeffZZ7NixAyqVc5crOTkZL7/8stPvQURE5JXUxt5ZNRWeHUcL4NEM1bPPPgtBEOr8SUtLs3hNdnY2Ro4cibvvvhvTp0+3e95t27ZhypQpWL58OXr0MDRKKy0txQMPPIDly5ejdevWdl9nTafT4b777sPLL7+Ma665xunPNXfuXBQXF8s/586x0ywREbVAPsaAqpYBVX2cqqG67rrrXDupIODbb79FmzZt6jzuwoULKCwsrPOYDh06yHfy5eTkYOjQobjhhhuQkpIChcI2Hty+fTtGjRqFRYsWYcaMGfL2Q4cOoW/fvlAqTYV1er2hk61CoUB6ejo6drQszisqKkJYWJjNa0RRhFKpxM8//4xbbrHtwmuNNVRERNQiXfwLeK8/oAkB5mZ5ejTNzu13+R06dAizZ89GYGBgvceKoojXX38d1dXV9R4bERGBiAjbpQvsyc7OxrBhw9CvXz+sWrXKbjCVmpqKpKQkLFiwwCKYAoCuXbviyJEjFtteeOEFlJaW4p133kF8vO2dEsHBwTav+eCDD7B161Z88cUXSEiwvyYXERHRFUHOULGGqj5O11DNmTMHkZGRTh371ltvNXhA9mRnZ2Po0KFo164dFi5ciAsXLsj7oqOjARim+ZKSkjBz5kyMGzcOeXl5AAC1Wo3w8HD4+vqiZ8+eFucNDQ0FAIvtc+fORXZ2Nj755BMoFAqb10RGRto9FxER0RVHbSxK12sBbQ2gcn6NwKuNUwHVmTNnnM4kAcDx48cRG2vb9bWhNm/ejIyMDGRkZCAuzrKnhjRjuXr1alRUVCA5ORnJycny/iFDhiA1NdXp98rNzUVW1tWX1iQiIrKhNlu+p6YMULmvN9eVpsX2oWpJWENFREQt1iutAX0t8MQxIMRxo9ArUbN0Sq+oqEBWVhZqamostvfu3buhpyQiIiJvo/YHqorZOqEeLgdUFy5cwJQpU/Djjz/a3a/T6Ro9KCIiIvISPgGGgIqF6XVyuQ/VrFmzUFRUhN9//x1+fn746aefsHr1anTu3BnffvttU4yRiIiIPIXNPZ3icoZq69at+Oabb9C/f38oFAq0a9cOt912G4KDg5GcnIxRo0Y1xTiJiIjIE9jc0ykuZ6jKy8vl9glhYWFyC4NevXrhjz/+cO/oiIiIyLOkO/1qOOVXF5cDqi5duiA9PR0AcO2112LZsmXIzs7Ghx9+iJiYGLcPkIiIiDyIGSqnuDzlN3PmTOTm5gIA5s2bh5EjR2LNmjVQq9VISUlx9/iIiIjIk5ihcorLAdXEiRPlx/369cPZs2eRlpaGtm3bOr3oMBEREbUQUkDFDFWdGtyHSuLv7+/y4slERETUQkhTfsxQ1cnpgOrJJ5906rhFixY1eDBERETkZdQMqJzhdEB18OBBi+c7d+5Ev3794OfnJ28TBMF9IyMiIiLP8+GUnzOcDqi2bdtm8TwoKAifffYZOnTo4PZBERERkZdgY0+nuNw2gYiIiK4ictsETvnVhQEVEREROSa3TWCGqi4MqIiIiMgxNvZ0itM1VIcPH7Z4Looi0tLSUFZWZrG9d+/e7hkZEREReR4bezrF6YCqT58+EAQBoijK25KSkgBA3i4IAnQ6nftHSURERJ7Bxp5OcTqgOnPmTFOOg4iIiLyRD+/yc4bTAVW7du2achxERETkjTjl5xSnitIPHz4MvV7v9EmPHTsGrVbb4EERERGRl2DbBKc4FVD17dsXhYWFTp80MTERWVlZDR4UEREReQmpsadeC2hrPDsWL+bUlJ8oinjxxRfh7+/v1ElranjBiYiIrgjS0jOAIUulUntuLF7MqYDq5ptvRnp6utMnTUxMtFjjj4iIiFoolRpQqAwZqpoKwC/M0yPySk4FVKmpqU08DCIiIvJaPgFAdTFbJ9SBndKJiIiobvICySxMd4QBFREREdWNy8/UiwEVERER1Y0LJNeLARURERHVTV5+hlN+jrgcUJWX82ISERFdVXxYQ1UflwOqqKgoTJ06FTt37myK8RAREZG3YVF6vVwOqD799FNcunQJt9xyC6655hq8/vrryMnJaYqxERERkTeQmnuyKN0hlwOqsWPHYsOGDcjOzsYjjzyCzz77DO3atUNSUhK++uorruFHRER0pZEzVAyoHGlwUXpERASefPJJHD58GIsWLcIvv/yCf/zjH4iNjcVLL72EigpedCIioisCF0iul1Od0u3Jz8/H6tWrkZKSgrNnz+If//gHpk2bhvPnz2PBggXYs2cPfv75Z3eOlYiIiDyBbRPq5XJA9dVXX2HVqlXYtGkTunfvjn/+85+YOHEiQkND5WMGDRqEbt26uXOcRERE5Cls7FkvlwOqKVOm4N5778Vvv/2G66+/3u4xsbGxeP755xs9OCIiIvICcoaKU36OuBxQ5ebmwt/fv85j/Pz8MG/evAYPioiIiLyImnf51cflgEqr1aKkpMRmuyAI0Gg0UKvVbhkYEREReQkf3uVXH5cDqtDQUAiC4HB/XFwcJk+ejHnz5kGh4Mo2RERELZ485Vfm2XF4MZcDqpSUFDz//POYPHkyBgwYAADYu3cvVq9ejRdeeAEXLlzAwoULodFo8Nxzz7l9wERERNTMWJReL5cDqtWrV+Ott97CPffcI28bPXo0evXqhWXLlmHLli1o27YtXn31VQZUREREVwI29qyXy3Nyu3btQt++fW229+3bF7t37wYA3HTTTcjKymr86IiIiMjz5KVneJefIy4HVPHx8VixYoXN9hUrViA+Ph4AUFhYiLCwsMaPzigzMxPTpk1DQkIC/Pz80LFjR8ybNw81NTXyMampqRgzZgxiYmIQEBCAPn36YM2aNQ7PuXbtWgiCgLFjx9b7/tXV1Xj++efRrl07aDQatG/fHitXrnTHRyMiIvJ+zFDVy+Upv4ULF+Luu+/Gjz/+KPeh2r9/P9LS0vDFF18AAPbt24fx48e7bZBpaWnQ6/VYtmwZOnXqhKNHj2L69OkoLy/HwoULARgyZ71798YzzzyDqKgobNy4EZMmTUJISAiSkpIszpeZmYmnnnoKgwcPdur977nnHuTn52PFihXo1KkTcnNzodfr3fb5iIiIvJpUQ6WvBXS1gNLHs+PxQoIoiqKrL8rMzMSyZcuQnp4OAOjSpQsefvhhtG/f3t3jc+jNN9/E0qVLcfr0aYfHjBo1ClFRURbZJJ1Oh5tvvhlTp07Fjh07UFRUhA0bNjg8x08//YR7770Xp0+fRnh4eIPGWlJSgpCQEBQXFyM4OLhB5yAiIvIYbTXwn0jD42fOAn6hHh1Oc3Hl+9ulDFVtbS1GjhyJDz/8EMnJyY0aZGMVFxfXG+AUFxfbLIHzyiuvIDIyEtOmTcOOHTvqfZ9vv/0W/fv3xxtvvIH//ve/CAgIwJ133on58+fDz8/P7muqq6tRXV0tP7fXt4uIiKjFUKoBhQrQaw13+l0lAZUrXAqofHx8cPjw4aYai9MyMjKwZMkSebrPnnXr1mHfvn1YtmyZvG3nzp1YsWIFDh065PR7nT59Gjt37oSvry++/vprXLx4Ef/85z9RWFiIVatW2X1NcnIyXn75Zaffg4iIyKsJgqEwvbqYdVQOuFyUPnHiRLtF6Q3x7LPPQhCEOn/S0tIsXpOdnY2RI0fi7rvvxvTp0+2ed9u2bZgyZQqWL1+OHj16AABKS0vxwAMPYPny5WjdurXTY9Tr9RAEAWvWrMGAAQNwxx13YNGiRVi9ejUqKyvtvmbu3LkoLi6Wf86dO+f0+xEREXklqTCdd/rZ1aClZ1auXIlffvkF/fr1Q0BAgMX+RYsWOX2u2bNnY/LkyXUe06FDB/lxTk4Ohg0bhkGDBuGjjz6ye/z27dsxevRoLF68GJMmTZK3nzp1CpmZmRg9erS8TSosV6lUSE9PR8eOHW3OFxMTgzZt2iAkJETe1q1bN4iiiPPnz6Nz5842r9FoNNBoNHV+LiIiohZFXn6GAZU9LgdUR48exXXXXQcAOHnypMW+upaksSciIgIRERFOHZudnY1hw4ahX79+WLVqld1lbVJTU5GUlIQFCxZgxowZFvu6du2KI0eOWGx74YUXUFpainfeeUdu+WDtxhtvxPr161FWVobAwEAAhs+tUCgQFxfn1NiJiIhaPLZOqJPLAdW2bduaYhx1ys7OxtChQ9GuXTssXLgQFy5ckPdFR0fL40pKSsLMmTMxbtw45OXlAQDUajXCw8Ph6+uLnj17Wpw3NDQUACy2z507F9nZ2fjkk08AAPfddx/mz5+PKVOm4OWXX8bFixcxZ84cTJ061WFROhER0RWHzT3r1ODVizMyMrBp0ya5jqgB3RectnnzZmRkZGDLli2Ii4tDTEyM/CNZvXo1KioqkJycbLH/rrvucum9cnNzLbq8BwYGYvPmzSgqKkL//v1x//33Y/To0Xj33Xfd9vmIiIi8HjNUdXK5D1VhYSHuuecebNu2DYIg4K+//kKHDh0wdepUhIWF4a233mqqsbZY7ENFREQt3tr7gbSNwKi3gOsf8vRomoUr398uZ6ieeOIJ+Pj4ICsrC/7+/vL28ePH46effnJ9tEREROT91MYpP2ao7HK5hurnn3/Gpk2bbAqyO3fujLNnz7ptYERERORFpICqlgGVPS5nqMrLyy0yU5JLly6xVQAREdGVim0T6uRyQDV48GD5DjjA0CpBr9fjjTfewLBhw9w6OCIiIvISzFDVyeUpvzfeeAO33nor9u/fj5qaGjz99NM4duwYLl26hN9++60pxkhERESe5sO7/OricoaqZ8+eOHnyJG666SaMGTMG5eXluOuuu3Dw4EG7ncaJiIjoCiAXpZd5dhxeyuUMFQCEhITg+eefd/dYiIiIyFtJGSpO+dnVoICqqKgIe/fuRUFBgbwensR8/TwiIiK6QrCxZ51cDqi+++473H///SgrK0NwcLDF+n2CIDCgIiIiuhJx6Zk6uVxDNXv2bEydOhVlZWUoKirC5cuX5Z9Lly41xRiJiIjI05ihqpPLAVV2djYef/xxu72oiIiI6ArFGqo6uRxQjRgxAvv372+KsRAREZG3UgcafrOxp10u11CNGjUKc+bMwfHjx9GrVy/4+PhY7L/zzjvdNjgiIiLyEmpmqOoiiKIouvIChcJxUksQBOh0ukYP6krjymrVREREXqniEvBGguHxi4WAskGNAloUV76/Xb4a1m0SiIiI6CogNfYEDHf6KUM8NxYv5HINFREREV2FlGpAUBoes47KhtMB1R133IHi4mL5+euvv46ioiL5eWFhIbp37+7WwREREZGXEASz5WdYR2XN6YBq06ZNqK6ulp+/9tprFn2ntFot0tPT3Ts6IiIi8h5y6wRmqKw5HVBZ1667WMtORERELR2bezrEGioiIiJyDpefccjpgEoQBIt1+6RtREREdJVghsohp9smiKKIyZMnQ6PRAACqqqrwyCOPICDAEK2a11cRERHRFUgqSmdzTxtOB1QPPvigxfOJEyfaHDNp0qTGj4iIiIi8k1SUzrYJNpwOqFatWtWU4yAiIiJvxwyVQyxKJyIiIuf4sIbKEQZURERE5By5sWeZZ8fhhRhQERERkXPkxp7MUFljQEVERETOYdsEhxhQERERkXPUgYbfNaWeHYcXYkBFREREzgmIMPwuK/DsOLwQAyoiIiJyTlC04XdpnmfH4YUYUBEREZFzAqMMv8vyAVH07Fi8DAMqIiIico6UodJWAVXFnh2Ll2FARURERM7x8QN8QwyPy/I9OxYvw4CKiIiInBfIOip7GFARERGR84LM6qhIxoCKiIiInMcMlV0MqIiIiMh5zFDZxYCKiIiInMcMlV0MqIiIiMh5UusEZqgsMKAiIiIi50nNPZmhssCAioiIiJzHDJVdDKiIiIjIeVKGqroEqCn37Fi8SIsIqDIzMzFt2jQkJCTAz88PHTt2xLx581BTUyMfk5qaijFjxiAmJgYBAQHo06cP1qxZ4/Cca9euhSAIGDt2bL3vv2bNGlx77bXw9/dHTEwMpk6disLCQnd8NCIiopZFEwT4+Bsec9pP1iICqrS0NOj1eixbtgzHjh3D4sWL8eGHH+K5556Tj9m1axd69+6NL7/8EocPH8aUKVMwadIkbNy40eZ8mZmZeOqppzB48OB63/u3337DpEmTMG3aNBw7dgzr16/H3r17MX36dLd+RiIiohZBECwXSSYAgCCKLXO56DfffBNLly7F6dOnHR4zatQoREVFYeXKlfI2nU6Hm2++GVOnTsWOHTtQVFSEDRs2ODzHwoULsXTpUpw6dUretmTJEixYsADnz5+3+5rq6mpUV1fLz0tKShAfH4/i4mIEBwe78CmJiIi80MqRQNZu4B+rgJ53eXo0TaakpAQhISFOfX+3iAyVPcXFxQgPD3f5mFdeeQWRkZGYNm2aU++TmJiIc+fO4YcffoAoisjPz8cXX3yBO+64w+FrkpOTERISIv/Ex8c79V5EREQtAjNUNlpkQJWRkYElS5bg4YcfdnjMunXrsG/fPkyZMkXetnPnTqxYsQLLly93+r1uvPFGrFmzBuPHj4darUZ0dDRCQkLw/vvvO3zN3LlzUVxcLP+cO3fO6fcjIiLyekFs7mnNowHVs88+C0EQ6vxJS0uzeE12djZGjhyJu+++22Ed07Zt2zBlyhQsX74cPXr0AACUlpbigQcewPLly9G6dWunx3j8+HHMnDkTL730Eg4cOICffvoJmZmZeOSRRxy+RqPRIDg42OKHiIjoisEMlQ2P1lBduHCh3rvlOnToALVaDQDIycnB0KFDccMNNyAlJQUKhW08uH37dowaNQqLFi3CjBkz5O2HDh1C3759oVQq5W16vR4AoFAokJ6ejo4dO9qc74EHHkBVVRXWr18vb9u5cycGDx6MnJwcxMTE1Ps5XZmDJSIi8nqHPgM2/B/QYRgwaYOnR9NkXPn+VjXTmOyKiIhARESEU8dmZ2dj2LBh6NevH1atWmU3mEpNTUVSUhIWLFhgEUwBQNeuXXHkyBGLbS+88AJKS0vxzjvvOKxzqqiogEpleZmkoKyF1vMTERE1DjNUNjwaUDkrOzsbQ4cORbt27bBw4UJcuHBB3hcdbZjH3bZtG5KSkjBz5kyMGzcOeXmGeV21Wo3w8HD4+vqiZ8+eFucNDQ0FAIvtc+fORXZ2Nj755BMAwOjRozF9+nQsXboUI0aMQG5uLmbNmoUBAwYgNja2KT82ERGRd2INlY0WEVBt3rwZGRkZyMjIQFxcnMU+KUu0evVqVFRUIDk5GcnJyfL+IUOGIDU11en3ys3NRVZWlvx88uTJKC0txXvvvYfZs2cjNDQUt9xyCxYsWNC4D0VERNRSBRoDqspLgLYGUKk9Ox4v0GL7ULUkrKEiIqIriigC8yMAfS0w6ygQemW2B7oq+lARERGRh7Bbug0GVEREROS6IGNAxToqAAyoiIiIqCGkOqoyBlQAAyoiIiJqCDlDxSk/gAEVERERNQQzVBYYUBEREZHrmKGywICKiIiIXMcMlQUGVEREROQ6uVs6M1QAAyoiIiJqCCmgKi8A9DrPjsULMKAiIiIi1wVEAIICEPVA+UVPj8bjGFARERGR6xRKQ1AFsI4KDKiIiIiooQJ5p5+EARURERE1TBDv9JMwoCIiIqKGCeR6fhIGVERERNQwcusEBlQMqIiIiKhhpAxVGWuoGFARERFRw0gZqpIcz47DCzCgIiIiooaJ7G74nX8UqK3y7Fg8jAEVERERNUx4B8O0n64GyD7g6dF4FAMqIiIiahhBANomGh5n7fLsWDyMARURERE1XLsbDb/PMqAiIiIiaph2xgzVub2ATuvZsXgQAyoiIiJquMjugG8IUFMG5B329Gg8hgEVERERNZxCCcTfYHh8FU/7MaAiIiKixmk3yPA7a7dnx+FBDKiIiIioccwL0/V6z47FQxhQERERUePEXAuo/IDKS8DFdE+PxiMYUBEREVHjqNRA/PWGx1dpHRUDKiIiImq8tsY6KgZURERERA3UziygEkXPjsUDGFARERFR48VdDyhUQGkOUHTW06NpdgyoiIiIqPHU/kBsX8Pjq3DaT+XpARAREdEVot0g4Pw+IGOLoYN6WQFQlg9EdAHiB3h6dE2KARURERG5R9tBwG/vAEe/MPxIlBpgdhrgH+65sTUxTvkRERGRe7S/CQjvAAgKIDAKiO4FaIIBXbUhc3UFY4aKiIiI3EMTCDz2ByDqDWv8AcDX/wf8+Rlwfj9wzQjPjq8JMUNFRERE7iMIpmAKAOL6GX5n7/fMeJoJAyoiIiJqOm36G35nH7ii1/ljQEVERERNJ6oHoPIFqoqBS6c8PZomw4CKiIiImo7SB4jpY3h8/sqd9mNARURERE0rzjjtdwXf6ceAioiIiJqWFFBdwYXpLSKgyszMxLRp05CQkAA/Pz907NgR8+bNQ01NjXxMamoqxowZg5iYGAQEBKBPnz5Ys2aNxXlSUlIgCILFj6+vb73vn5qaiuuuuw4ajQadOnVCSkqKuz8iERHRlUsqTM8/BtRWenYsTaRF9KFKS0uDXq/HsmXL0KlTJxw9ehTTp09HeXk5Fi5cCADYtWsXevfujWeeeQZRUVHYuHEjJk2ahJCQECQlJcnnCg4ORnp6uvxcEIQ63/vMmTMYNWoUHnnkEaxZswZbtmzBQw89hJiYGIwYceX20yAiInKbkDhDo8+yfCD3T6DtDZ4ekdsJoiiKnh5EQ7z55ptYunQpTp8+7fCYUaNGISoqCitXrgRgyFDNmjULRUVFTr/PM888g++//x5Hjx6Vt917770oKirCTz/95NQ5SkpKEBISguLiYgQHBzv93kRERFeMz+8D0r8H/vYqMOhfnh6NU1z5/m4RU372FBcXIzy87jWB7B1TVlaGdu3aIT4+HmPGjMGxY8fqPMfu3bsxfPhwi20jRozA7t27Hb6muroaJSUlFj9ERERXtSu8wWeLDKgyMjKwZMkSPPzwww6PWbduHfbt24cpU6bI27p06YKVK1fim2++waeffgq9Xo9Bgwbh/PnzDs+Tl5eHqKgoi21RUVEoKSlBZaX9eeDk5GSEhITIP/Hx8S5+QiIioitM3PWG31do6wSPBlTPPvusTZG49U9aWprFa7KzszFy5EjcfffdmD59ut3zbtu2DVOmTMHy5cvRo0cPeXtiYiImTZqEPn36YMiQIfjqq68QERGBZcuWufVzzZ07F8XFxfLPuXPn3Hp+IiKiFie2LwABKD4HlOZ7ejRu59Gi9NmzZ2Py5Ml1HtOhQwf5cU5ODoYNG4ZBgwbho48+snv89u3bMXr0aCxevBiTJk2q89w+Pj7o27cvMjIyHB4THR2N/HzLP/j8/HwEBwfDz8/P7ms0Gg00Gk2d701ERHRV0QQBkd2AguOGab+uozw9IrfyaEAVERGBiIgIp47Nzs7GsGHD0K9fP6xatQoKhW1yLTU1FUlJSViwYAFmzJhR7zl1Oh2OHDmCO+64w+ExiYmJ+OGHHyy2bd68GYmJiU6Nm4iIiIza9DMEVOevvICqRdRQZWdnY+jQoWjbti0WLlyICxcuIC8vD3l5efIx27Ztw6hRo/D4449j3Lhx8v5Lly7Jx7zyyiv4+eefcfr0afzxxx+YOHEizp49i4ceekg+Zu7cuRaZrUceeQSnT5/G008/jbS0NHzwwQdYt24dnnjiieb58ERERFeKK7jBZ4voQ7V582ZkZGQgIyMDcXFxFvukrg+rV69GRUUFkpOTkZycLO8fMmQIUlNTAQCXL1/G9OnTkZeXh7CwMPTr1w+7du1C9+7d5eNzc3ORlZUlP09ISMD333+PJ554Au+88w7i4uLw8ccfswcVERGRq6QGn9kHAb0OUCg9Ox43arF9qFoS9qEiIiKCIYhKjgdqy4F/7jHUVHmxq6IPFREREbUwCqUpiLp40rNjcTMGVERERNR8WnU0/C50fId9S8SAioiIiJpPq06G34WOl45riRhQERERUfMJN/aXvHTKs+NwMwZURERE1Hw45UdERETUSOHGgKr8AlBV4tmxuBEDKiIiImo+vsFAQKTh8RU07ceAioiIiJqXPO3HgIqIiIioYRhQERERETWSVEfFKT8iIiKiBroC7/RjQEVERETNS27uyQwVERERUcOEJRh+VxUBFZc8OhR3YUBFREREzUvtDwS3MTy+Qqb9GFARERFR87vC7vRjQEVERETN7wq7048BFRERETW/K+xOPwZURERE1PyusDv9GFARERFR85On/E4DoujZsbgBAyoiIiJqfmHtAUEB1JQBZfmeHk2jMaAiIiKi5qdSA6FtDY+vgGk/BlRERETkGVfQnX4MqIiIiMgz5ML0ln+nHwMqIiIi8owrqLknAyoiIiLyDPM7/VxReAqoKXf/eBqBARURERF5RiuzgEqvd+41Wb8DS/oBG/7ZdONqAAZURERE5Bkh8YDCB9BWASXZzr3mxLcARCDte6C6tEmH5woGVEREROQZSpWhHxUAFP7l3Gsydxh+62uBU9uaZFgNwYCKiIiIPCe2r+H3vhX1H1txCcg9bHp+8qemGVMDMKAiIiIizxk8GxCUQNpGIHNn3cee3QVABJRqw/OTm5yvvWpiDKiIiIjIcyK7Av0mGx5veq7uAOnMr4bf104ANMFAxUUg548mH6IzGFARERGRZw2dC6iDgNw/gSPrHB8n1U91vAXodKvhcfqPTT8+JzCgIiIiIs8KjABunm14/MvLQE2F7TFlF4CC44bH7QcD19xueHxyU/OMsR4MqIiIiMjzBv4fENIWKM0Bdr9vu1/KTkX1BAJaAZ2GA4ICyD8CFJ9v3rHawYCKiIiIPM/HFxg+z/B452KgNM9yv1Q/lXCz4XdAKyBugOGxF9ztx4CKiIiIvEPPcUCb/kBtuaFA3ZyUoWo/2LSty0jDby+Y9mNARURERN5BEIA73jRM5R390lRwXpwNFGYYtrcbZDr+GmNAdXq7x9f2Y0BFRERE3qPNdUDivwyPNz4BVBWbslMxfQC/UNOxEV2B0LaArtoQVHkQAyoiIiLyLsOeA8I7AKW5wM8vAmeMAVXCYMvjBMHsbj/P1lExoCIiIiLv4uMH3Pme4fEfq4Hj3xgeSwXp5q4ZYfjt4a7pDKiIiIjI+7S/Eeg/zfC4phRQqID4G+wcdxMQex1w7b2Atqp5x2hG5bF3JiIiIqrL8H8bMk8l5w13/2kCbY9RaYAZ25p9aNZaRIYqMzMT06ZNQ0JCAvz8/NCxY0fMmzcPNTU18jGpqakYM2YMYmJiEBAQgD59+mDNmjUW50lJSYEgCBY/vr6+db73V199hdtuuw0REREIDg5GYmIiNm3y/O2ZREREVzzfYODvHwJhCcCA6Z4eTZ1aRIYqLS0Ner0ey5YtQ6dOnXD06FFMnz4d5eXlWLhwIQBg165d6N27N5555hlERUVh48aNmDRpEkJCQpCUlCSfKzg4GOnp6fJzQRDqfO9ff/0Vt912G1577TWEhoZi1apVGD16NH7//Xf07du3aT4wERERGSQMBmYe8vQo6iWIoih6ehAN8eabb2Lp0qU4ffq0w2NGjRqFqKgorFy5EoAhQzVr1iwUFRU16r179OiB8ePH46WXXnLq+JKSEoSEhKC4uBjBwcGNem8iIiJqHq58f7eIKT97iouLER4e7vIxZWVlaNeuHeLj4zFmzBgcO3bMpffV6/UoLS2t872rq6tRUlJi8UNERERXrhYZUGVkZGDJkiV4+OGHHR6zbt067Nu3D1OmTJG3denSBStXrsQ333yDTz/9FHq9HoMGDcL5884vqrhw4UKUlZXhnnvucXhMcnIyQkJC5J/4+Hinz09EREQtj0en/J599lksWLCgzmNOnDiBrl27ys+zs7MxZMgQDB06FB9//LHd12zbtg1JSUlYunQpJk2a5PDctbW16NatGyZMmID58+fXO97PPvsM06dPxzfffIPhw4c7PK66uhrV1dXy85KSEsTHx3PKj4iIqAVxZcrPo0Xps2fPxuTJk+s8pkOHDvLjnJwcDBs2DIMGDcJHH31k9/jt27dj9OjRWLx4cZ3BFAD4+Pigb9++yMjIqHesa9euxUMPPYT169fXGUwBgEajgUajqfecREREdGXwaEAVERGBiIgIp47Nzs7GsGHD0K9fP6xatQoKhe1sZWpqKpKSkrBgwQLMmDGj3nPqdDocOXIEd9xxR53Hff7555g6dSrWrl2LUaNGOTVeIiIiunq0iLYJ2dnZGDp0KNq1a4eFCxfiwoUL8r7o6GgApmm+mTNnYty4ccjLywMAqNVquYD8lVdewQ033IBOnTqhqKgIb775Js6ePYuHHnpIPt/cuXORnZ2NTz75BIBhmu/BBx/EO++8g4EDB8rn9fPzQ0hISLN8fiIiIvJuLaIoffPmzcjIyMCWLVsQFxeHmJgY+UeyevVqVFRUIDk52WL/XXfdJR9z+fJlTJ8+Hd26dcMdd9yBkpIS7Nq1C927d5ePyc3NRVZWlvz8o48+glarxaOPPmpx3pkzZzbPhyciIiKv12L7ULUk7ENFRETU8lwVfaiIiIiIvAUDKiIiIqJGYkBFRERE1EgMqIiIiIgaiQEVERERUSO1iD5ULZ10IyUXSSYiImo5pO9tZxoiMKBqBqWlpQDARZKJiIhaoNLS0nqbebMPVTPQ6/XIyclBUFAQBEFw67mlhZfPnTvHHleNxGvpHryO7sNr6T68lu5zNV1LURRRWlqK2NhYu0vemWOGqhkoFArExcU16XsEBwdf8X+xmwuvpXvwOroPr6X78Fq6z9VyLZ1dZo5F6URERESNxICKiIiIqJEYULVwGo0G8+bNg0aj8fRQWjxeS/fgdXQfXkv34bV0H15L+1iUTkRERNRIzFARERERNRIDKiIiIqJGYkBFRERE1EgMqIiIiIgaiQFVC/b++++jffv28PX1xcCBA7F3715PD8nrJScn4/rrr0dQUBAiIyMxduxYpKenWxxTVVWFRx99FK1atUJgYCDGjRuH/Px8D424ZXj99dchCAJmzZolb+N1dF52djYmTpyIVq1awc/PD7169cL+/fvl/aIo4qWXXkJMTAz8/PwwfPhw/PXXXx4csXfS6XR48cUXkZCQAD8/P3Ts2BHz58+3WIeN19K+X3/9FaNHj0ZsbCwEQcCGDRss9jtz3S5duoT7778fwcHBCA0NxbRp01BWVtaMn8KzGFC1UP/73//w5JNPYt68efjjjz9w7bXXYsSIESgoKPD00Lza9u3b8eijj2LPnj3YvHkzamtr8be//Q3l5eXyMU888QS+++47rF+/Htu3b0dOTg7uuusuD47au+3btw/Lli1D7969LbbzOjrn8uXLuPHGG+Hj44Mff/wRx48fx1tvvYWwsDD5mDfeeAPvvvsuPvzwQ/z+++8ICAjAiBEjUFVV5cGRe58FCxZg6dKleO+993DixAksWLAAb7zxBpYsWSIfw2tpX3l5Oa699lq8//77dvc7c93uv/9+HDt2DJs3b8bGjRvx66+/YsaMGc31ETxPpBZpwIAB4qOPPio/1+l0YmxsrJicnOzBUbU8BQUFIgBx+/btoiiKYlFRkejj4yOuX79ePubEiRMiAHH37t2eGqbXKi0tFTt37ixu3rxZHDJkiDhz5kxRFHkdXfHMM8+IN910k8P9er1ejI6OFt988015W1FRkajRaMTPP/+8OYbYYowaNUqcOnWqxba77rpLvP/++0VR5LV0FgDx66+/lp87c92OHz8uAhD37dsnH/Pjjz+KgiCI2dnZzTZ2T2KGqgWqqanBgQMHMHz4cHmbQqHA8OHDsXv3bg+OrOUpLi4GAISHhwMADhw4gNraWotr27VrV7Rt25bX1o5HH30Uo0aNsrheAK+jK7799lv0798fd999NyIjI9G3b18sX75c3n/mzBnk5eVZXMuQkBAMHDiQ19LKoEGDsGXLFpw8eRIA8Oeff2Lnzp24/fbbAfBaNpQz12337t0IDQ1F//795WOGDx8OhUKB33//vdnH7AlcHLkFunjxInQ6HaKioiy2R0VFIS0tzUOjann0ej1mzZqFG2+8ET179gQA5OXlQa1WIzQ01OLYqKgo5OXleWCU3mvt2rX4448/sG/fPpt9vI7OO336NJYuXYonn3wSzz33HPbt24fHH38carUaDz74oHy97P33zmtp6dlnn0VJSQm6du0KpVIJnU6HV199Fffffz8A8Fo2kDPXLS8vD5GRkRb7VSoVwsPDr5pry4CKrlqPPvoojh49ip07d3p6KC3OuXPnMHPmTGzevBm+vr6eHk6Lptfr0b9/f7z22msAgL59++Lo0aP48MMP8eCDD3p4dC3LunXrsGbNGnz22Wfo0aMHDh06hFmzZiE2NpbXkpocp/xaoNatW0OpVNrcMZWfn4/o6GgPjapl+de//oWNGzdi27ZtiIuLk7dHR0ejpqYGRUVFFsfz2lo6cOAACgoKcN1110GlUkGlUmH79u149913oVKpEBUVxevopJiYGHTv3t1iW7du3ZCVlQUA8vXif+/1mzNnDp599lnce++96NWrFx544AE88cQTSE5OBsBr2VDOXLfo6Gibm6K0Wi0uXbp01VxbBlQtkFqtRr9+/bBlyxZ5m16vx5YtW5CYmOjBkXk/URTxr3/9C19//TW2bt2KhIQEi/39+vWDj4+PxbVNT09HVlYWr62ZW2+9FUeOHMGhQ4fkn/79++P++++XH/M6OufGG2+0ad1x8uRJtGvXDgCQkJCA6Ohoi2tZUlKC33//ndfSSkVFBRQKy681pVIJvV4PgNeyoZy5bomJiSgqKsKBAwfkY7Zu3Qq9Xo+BAwc2+5g9wtNV8dQwa9euFTUajZiSkiIeP35cnDFjhhgaGirm5eV5emhe7f/+7//EkJAQMTU1VczNzZV/Kioq5GMeeeQRsW3btuLWrVvF/fv3i4mJiWJiYqIHR90ymN/lJ4q8js7au3evqFKpxFdffVX866+/xDVr1oj+/v7ip59+Kh/z+uuvi6GhoeI333wjHj58WBwzZoyYkJAgVlZWenDk3ufBBx8U27RpI27cuFE8c+aM+NVXX4mtW7cWn376afkYXkv7SktLxYMHD4oHDx4UAYiLFi0SDx48KJ49e1YUReeu28iRI8W+ffuKv//+u7hz506xc+fO4oQJEzz1kZodA6oWbMmSJWLbtm1FtVotDhgwQNyzZ4+nh+T1ANj9WbVqlXxMZWWl+M9//lMMCwsT/f39xb///e9ibm6u5wbdQlgHVLyOzvvuu+/Enj17ihqNRuzatav40UcfWezX6/Xiiy++KEZFRYkajUa89dZbxfT0dA+N1nuVlJSIM2fOFNu2bSv6+vqKHTp0EJ9//nmxurpaPobX0r5t27bZ/X/jgw8+KIqic9etsLBQnDBhghgYGCgGBweLU6ZMEUtLSz3waTxDEEWzFrJERERE5DLWUBERERE1EgMqIiIiokZiQEVERETUSAyoiIiIiBqJARURERFRIzGgIiIiImokBlREREREjcSAioiIiKiRGFARUYvw73//G3369GnUOTIzMyEIAg4dOuSWMTkydOhQzJo1q0nfg4i8CwMqInKLc+fOYerUqYiNjYVarUa7du0wc+ZMFBYWunwuQRCwYcMGi21PPfWUxeKsDREfH4/c3Fz07NmzUeeRpKamQhAEFBUVWWz/6quvMH/+fLe8R0M0V+BIRCYMqIio0U6fPo3+/fvjr7/+wueff46MjAx8+OGH2LJlCxITE3Hp0qVGv0dgYCBatWrVqHMolUpER0dDpVI1ejx1CQ8PR1BQUJO+BxF5FwZURNRojz76KNRqNX7++WcMGTIEbdu2xe23345ffvkF2dnZeP755+Vj27dvj/nz52PChAkICAhAmzZt8P7771vsB4C///3vEARBfm495Td58mSMHTsWr732GqKiohAaGopXXnkFWq0Wc+bMQXh4OOLi4rBq1Sr5NdaZm8mTJ0MQBJuf1NRUAMB///tf9O/fH0FBQYiOjsZ9992HgoIC+VzDhg0DAISFhUEQBEyePBmA7ZTf5cuXMWnSJISFhcHf3x+33347/vrrL3l/SkoKQkNDsWnTJnTr1g2BgYEYOXIkcnNzHV7zy5cv4/7770dERAT8/PzQuXNn+bMmJCQAAPr27QtBEDB06FD5dR9//DG6desGX19fdO3aFR988IHN9Vm7di0GDRoEX19f9OzZE9u3b3fqfYmuap5enZmIWrbCwkJREATxtddes7t/+vTpYlhYmKjX60VRFMV27dqJQUFBYnJyspieni6+++67olKpFH/++WdRFEWxoKBABCCuWrVKzM3NFQsKCkRRFMV58+aJ1157rXzeBx98UAwKChIfffRRMS0tTVyxYoUIQBwxYoT46quviidPnhTnz58v+vj4iOfOnRNFURTPnDkjAhAPHjwoiqIoFhUVibm5ufLPzJkzxcjISDE3N1cURVFcsWKF+MMPP4inTp0Sd+/eLSYmJoq33367KIqiqNVqxS+//FIEIKanp4u5ubliUVGRKIqiOGTIEHHmzJnyWO+8806xW7du4q+//ioeOnRIHDFihNipUyexpqZGFEVRXLVqlejj4yMOHz5c3Ldvn3jgwAGxW7du4n333efwuj/66KNinz59xH379olnzpwRN2/eLH777beiKIri3r17RQDiL7/8Iubm5oqFhYWiKIrip59+KsbExIhffvmlePr0afHLL78Uw8PDxZSUFIvrExcXJ37xxRfi8ePHxYceekgMCgoSL168WO/7El3NGFARUaPs2bNHBCB+/fXXdvcvWrRIBCDm5+eLomgIqEaOHGlxzPjx4+VARRRFu+ezF1C1a9dO1Ol08rYuXbqIgwcPlp9rtVoxICBA/Pzzz0VRtA2ozH355Zeir6+vuHPnToefdd++fSIAsbS0VBRFUdy2bZsIQLx8+bLFceYB1cmTJ0UA4m+//Sbvv3jxoujn5yeuW7dOFEVDQAVAzMjIkI95//33xaioKIdjGT16tDhlyhS7+xx9zo4dO4qfffaZxbb58+eLiYmJFq97/fXX5f21tbViXFycuGDBgnrfl+hqxik/InILURSdPjYxMdHm+YkTJ1x+zx49ekChMP1vLCoqCr169ZKfK5VKtGrVSp6mc+TgwYN44IEH8N577+HGG2+Utx84cACjR49G27ZtERQUhCFDhgAAsrKynB7jiRMnoFKpMHDgQHlbq1at0KVLF4vP7O/vj44dO8rPY2Ji6hz3//3f/2Ht2rXo06cPnn76aezatavOcZSXl+PUqVOYNm0aAgMD5Z///Oc/OHXqlMWx5n8+KpUK/fv3l8fq6vsSXS0YUBFRo3Tq1AmCIDgMiE6cOIGwsDBERES4/b19fHwsnguCYHebXq93eI68vDzceeedeOihhzBt2jR5e3l5OUaMGIHg4GCsWbMG+/btw9dffw0AqKmpceOnMLA37rqC1Ntvvx1nz57FE088gZycHNx666146qmnHB5fVlYGAFi+fDkOHTok/xw9ehR79uxxepyuvi/R1YIBFRE1SqtWrXDbbbfhgw8+QGVlpcW+vLw8rFmzBuPHj4cgCPJ26y/wPXv2oFu3bvJzHx8f6HS6ph04gKqqKowZMwZdu3bFokWLLPalpaWhsLAQr7/+OgYPHoyuXbvaZIzUajUA1DnWbt26QavV4vfff5e3FRYWIj09Hd27d2/U+CMiIvDggw/i008/xdtvv42PPvrI4biioqIQGxuL06dPo1OnThY/UhG7xPzPR6vV4sCBAxZ/Po7el+hq1rT3DhPRVeG9997DoEGDMGLECPznP/9BQkICjh07hjlz5qBNmzZ49dVXLY7/7bff8MYbb2Ds2LHYvHkz1q9fj++//17e3759e2zZsgU33ngjNBoNwsLCmmTcDz/8MM6dO4ctW7bgwoUL8vbw8HC0bdsWarUaS5YswSOPPIKjR4/a9JZq164dBEHAxo0bcccdd8DPzw+BgYEWx3Tu3BljxozB9OnTsWzZMgQFBeHZZ59FmzZtMGbMmAaP/aWXXkK/fv3Qo0cPVFdXY+PGjXLQExkZCT8/P/z000+Ii4uDr68vQkJC8PLLL+Pxxx9HSEgIRo4cierqauzfvx+XL1/Gk08+KZ/7/fffR+fOndGtWzcsXrwYly9fxtSpU+t9X6KrGTNURNRonTt3xv79+9GhQwfcc8896NixI2bMmIFhw4Zh9+7dCA8Ptzh+9uzZ2L9/P/r27Yv//Oc/WLRoEUaMGCHvf+utt7B582bEx8ejb9++TTbu7du3Izc3F927d0dMTIz8s2vXLkRERCAlJQXr169H9+7d8frrr2PhwoUWr2/Tpg1efvllPPvss4iKisK//vUvu++zatUq9OvXD0lJSUhMTIQoivjhhx9spvlcoVarMXfuXPTu3Rs333wzlEol1q5dC8BQ9/Tuu+9i2bJliI2NlQO3hx56CB9//DFWrVqFXr16YciQIUhJSbHJUL3++ut4/fXXce2112Lnzp349ttv0bp163rfl+hqJoiuVJISETVS+/btMWvWLC7N4oUyMzORkJCAgwcPNnqZH6KrDTNURERERI3EgIqIiIiokTjlR0RERNRIzFARERERNRIDKiIiIqJGYkBFRERE1EgMqIiIiIgaiQEVERERUSMxoCIiIiJqJAZURERERI3EgIqIiIiokf4fRILvMqqMS0kAAAAASUVORK5CYII=", "text/plain": [ "
" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "# Plot the energies.\n", "\n", "vqe_y = vqe_energies\n", "vqe_x = list(range(len(vqe_y)))\n", "plt.plot(vqe_x, vqe_y, label=\"VQE\")\n", "\n", "afqmc_y = list(qmc_data[\"ETotal\"])\n", "afqmc_x = [i + vqe_x[-1] for i in list(range(len(afqmc_y)))]\n", "plt.plot(afqmc_x, afqmc_y, label=\"AFQMC\")\n", "\n", "plt.xlabel(\"Optimization steps\")\n", "plt.ylabel(\"Energy [Ha]\")\n", "plt.legend()" ] }, { "cell_type": "code", "execution_count": 11, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "CUDA-Q Version latest (https://github.com/NVIDIA/cuda-quantum 176f1e7df8a58c2dc3d6b1b47bf7f63b4b8d3b63)\n" ] } ], "source": [ "print(cudaq.__version__)" ] } ], "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.10.12" } }, "nbformat": 4, "nbformat_minor": 2 }