{ "cells": [ { "cell_type": "markdown", "source": [ "# Tutorial" ], "metadata": {} }, { "cell_type": "markdown", "source": [ "DFTK is a Julia package for playing with plane-wave\n", "density-functional theory algorithms. In its basic formulation it\n", "solves periodic Kohn-Sham equations.\n", "\n", "This document provides an overview of the structure of the code\n", "and how to access basic information about calculations.\n", "Basic familiarity with the concepts of plane-wave density functional theory\n", "is assumed throughout. Feel free to take a look at the\n", "[Periodic problems](https://docs.dftk.org/stable/guide/periodic_problems/)\n", "or the\n", "[Introductory resources](https://docs.dftk.org/stable/guide/introductory_resources/)\n", "chapters for some introductory material on the topic.\n", "\n", "!!! note \"Convergence parameters in the documentation\"\n", " We use rough parameters in order to be able\n", " to automatically generate this documentation very quickly.\n", " Therefore results are far from converged.\n", " Tighter thresholds and larger grids should be used for more realistic results.\n", "\n", "For our discussion we will use the classic example of\n", "computing the LDA ground state of the\n", "[silicon crystal](https://www.materialsproject.org/materials/mp-149).\n", "Performing such a calculation roughly proceeds in three steps." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "using DFTK\n", "using Plots\n", "using Unitful\n", "using UnitfulAtomic\n", "\n", "# 1. Define lattice and atomic positions\n", "a = 5.431u\"angstrom\" # Silicon lattice constant\n", "lattice = a / 2 * [[0 1 1.]; # Silicon lattice vectors\n", " [1 0 1.]; # specified column by column\n", " [1 1 0.]];" ], "metadata": {}, "execution_count": 1 }, { "cell_type": "markdown", "source": [ "By default, all numbers passed as arguments are assumed to be in atomic\n", "units. Quantities such as temperature, energy cutoffs, lattice vectors, and\n", "the k-point grid spacing can optionally be annotated with Unitful units,\n", "which are automatically converted to the atomic units used internally. For\n", "more details, see the [Unitful package\n", "documentation](https://painterqubits.github.io/Unitful.jl/stable/) and the\n", "[UnitfulAtomic.jl package](https://github.com/sostock/UnitfulAtomic.jl)." ], "metadata": {} }, { "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "n Energy log10(ΔE) log10(Δρ) Diag\n", "--- --------------- --------- --------- ----\n", " 1 -7.900411135525 -0.70 4.1\n", " 2 -7.905004191573 -2.34 -1.52 1.0\n", " 3 -7.905210674200 -3.69 -2.52 3.0\n", " 4 -7.905211513721 -6.08 -3.33 2.8\n", " 5 -7.905211530290 -7.78 -4.57 1.4\n", " 6 -7.905211531399 -8.95 -5.22 3.6\n" ] } ], "cell_type": "code", "source": [ "# Load HGH pseudopotential for Silicon\n", "Si = ElementPsp(:Si, psp=load_psp(\"hgh/lda/Si-q4\"))\n", "\n", "# Specify type and positions of atoms\n", "atoms = [Si, Si]\n", "positions = [ones(3)/8, -ones(3)/8]\n", "\n", "# 2. Select model and basis\n", "model = model_LDA(lattice, atoms, positions)\n", "kgrid = [4, 4, 4] # k-point grid (Regular Monkhorst-Pack grid)\n", "Ecut = 7 # kinetic energy cutoff\n", "# Ecut = 190.5u\"eV\" # Could also use eV or other energy-compatible units\n", "basis = PlaneWaveBasis(model; Ecut, kgrid)\n", "# Note the implicit passing of keyword arguments here:\n", "# this is equivalent to PlaneWaveBasis(model; Ecut=Ecut, kgrid=kgrid)\n", "\n", "# 3. Run the SCF procedure to obtain the ground state\n", "scfres = self_consistent_field(basis, tol=1e-8);" ], "metadata": {}, "execution_count": 2 }, { "cell_type": "markdown", "source": [ "That's it! Now you can get various quantities from the result of the SCF.\n", "For instance, the different components of the energy:" ], "metadata": {} }, { "outputs": [ { "output_type": "execute_result", "data": { "text/plain": "Energy breakdown (in Ha):\n Kinetic 3.1020988 \n AtomicLocal -2.1987894\n AtomicNonlocal 1.7296109 \n Ewald -8.3979253\n PspCorrection -0.2946254\n Hartree 0.5530407 \n Xc -2.3986218\n\n total -7.905211531399" }, "metadata": {}, "execution_count": 3 } ], "cell_type": "code", "source": [ "scfres.energies" ], "metadata": {}, "execution_count": 3 }, { "cell_type": "markdown", "source": [ "Eigenvalues:" ], "metadata": {} }, { "outputs": [ { "output_type": "execute_result", "data": { "text/plain": "7×8 Matrix{Float64}:\n -0.176942 -0.147441 -0.0911694 … -0.101219 -0.0239772 -0.0184082\n 0.261073 0.116914 0.004825 0.0611641 -0.0239772 -0.0184082\n 0.261073 0.232989 0.216733 0.121636 0.155531 0.117747\n 0.261073 0.232989 0.216733 0.212134 0.155531 0.117747\n 0.354532 0.335109 0.317102 0.350436 0.285692 0.417258\n 0.354532 0.389828 0.384601 … 0.436926 0.285692 0.417508\n 0.354532 0.389828 0.384601 0.449273 0.627568 0.443806" }, "metadata": {}, "execution_count": 4 } ], "cell_type": "code", "source": [ "hcat(scfres.eigenvalues...)" ], "metadata": {}, "execution_count": 4 }, { "cell_type": "markdown", "source": [ "`eigenvalues` is an array (indexed by k-points) of arrays (indexed by\n", "eigenvalue number). The \"splatting\" operation `...` calls `hcat`\n", "with all the inner arrays as arguments, which collects them into a\n", "matrix.\n", "\n", "The resulting matrix is 7 (number of computed eigenvalues) by 8\n", "(number of irreducible k-points). There are 7 eigenvalues per\n", "k-point because there are 4 occupied states in the system (4 valence\n", "electrons per silicon atom, two atoms per unit cell, and paired\n", "spins), and the eigensolver gives itself some breathing room by\n", "computing some extra states (see `n_ep_extra` argument to\n", "`self_consistent_field`). There are only 8 k-points (instead of\n", "4x4x4) because symmetry has been used to reduce the amount of\n", "computations to just the irreducible k-points (see\n", "[Crystal symmetries](https://docs.dftk.org/stable/developer/symmetries/)\n", "for details).\n", "\n", "We can check the occupations ..." ], "metadata": {} }, { "outputs": [ { "output_type": "execute_result", "data": { "text/plain": "7×8 Matrix{Float64}:\n 2.0 2.0 2.0 2.0 2.0 2.0 2.0 2.0\n 2.0 2.0 2.0 2.0 2.0 2.0 2.0 2.0\n 2.0 2.0 2.0 2.0 2.0 2.0 2.0 2.0\n 2.0 2.0 2.0 2.0 2.0 2.0 2.0 2.0\n 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0\n 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0\n 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0" }, "metadata": {}, "execution_count": 5 } ], "cell_type": "code", "source": [ "hcat(scfres.occupation...)" ], "metadata": {}, "execution_count": 5 }, { "cell_type": "markdown", "source": [ "... and density, where we use that the density objects in DFTK are\n", "indexed as ρ[iσ, ix, iy, iz], i.e. first in the spin component and then\n", "in the 3-dimensional real-space grid." ], "metadata": {} }, { "outputs": [ { "output_type": "execute_result", "data": { "text/plain": "Plot{Plots.GRBackend() n=1}", "image/png": "iVBORw0KGgoAAAANSUhEUgAAAlgAAAGQCAIAAAD9V4nPAAAABmJLR0QA/wD/AP+gvaeTAAAgAElEQVR4nOzdd2BTVdsA8OfcJN17pS3QXYZQOmiBFhDKkKWAoigCovK+OBAERUFli4ofr0y3ryj4KksQcaFAy2wp0DJKB1sKpE26m+4k93x/RGspLZQ2N3c9v7+a9DQ5Oc05zz3zEkopIIQQQnLF8J0BhBBCiE8YCBFCCMkaBkKEEEKyhoEQIYSQrGEgRAghJGsYCBFCCMkaBkKEEEKyhoEQIYSQrGEgRAghJGsYCBFCCMma1ALh9evXdTpdKxObTCZOMyNJWGhtgIXWBlhobYCF1jZSC4QrV67cunVrKxNXV1dzmhlJwkJrAyy0NsBCawMstLaRWiBECCGE7gkGQoQQQrKGgRAhhJCsYSBECCEkaxgIEUIIyRoGQoQQQrIm30D4v81bH//XS/+3Zr3RaOQ7LwjxICsr64l/zZg++3WNRsN3XhDik0wD4c5dP760fvtvXWctTbrx5rJ3+c4OQtZWXl4+7NEpWz3G/9cUn/jQY3xnByE+yTQQ7t57sHzADOjUs3rY6zt+P8B3dhCytpNnsysCE6DzABr5oJZ1LCgq4TtHCPFGpoEwMT7W+cQ3UJavPPJFsV9c7C7jf8+zVThEimTgbAl98ajpsXPBpkvHQJMDl1Pr9CVRe5zeOGG6qqd85w4hHsg0EE6d/OSicb16/frCzNCa/C1LV8Qp/rhBAzYbnjtiOluCbQGSoDoTbL/KDvvNOHKPycUGzjyl3v/N+oEnlz/451dnf9t8ZIyKAvTdbRz2m3H7VdbI8p1dhKyIUCqpdn/WrFnh4eEzZ85sTWK9Xu/s7NzwML8aNl1kP85hve1geldmShhjr+Qso6LVpNBQa/BbaJcr6Bfn2a8vsBEeZHpX5uFARtnCBXCdCXbnsZ/nstmlMCWczLiP6eRIrJvZf+A3rQ2w0NpGpj3CZvk5wLxI5urjyhVxip/y2JCthvk4WIREi6Ww7yadsN/U50djjRFSxij3jlQ+FtxiFAQAWwU8FszsHancN0pRY4SIHcaH/jDuuymti2WEboNdnqYYAkM7kKEdlJcq6H/Ps313G3ve7ToaIUEpqIGNF9hPclgvO5jeldk4UHWvYxvd3MjaeMXyWMXmy+zcNFONCZ7tzPyrK+Npy02OEeIVNu0tCnMhK+IUeU+opndlPs9lg7Ya558wXa/Ci2MkXEcK6IT9pu7fG67o6Y8PKE6OU07v2vYRfmcVTO/KnH5E+V2i4oqehm8zTNhv2ncTqwCSGuwR3oV5sOixYCanjH6aw0btNA7xZ6Z3ZYZ0ILxNniB0q/J62HqFXZfFshSe68psuF/lpLLk6/fyIp/1V7zfW7HxAvv8UZObDUzvykwKYxyx/UCSgD3C1jIPFv35hGpoBzI3zdR1u/H9M2xx3V+/NZlMJSW4Ewtxrr6+vry8vOFhehF97ogpeKth3026pq8i+1Hlyz0Yy0bBBm428HIP5sJjyhVxin03/1plndlolXVJSQneIR1ZQWFhoWWXeeIV3b0xDxZN78qkF9HPc9nwbYah/kxCRer/vTHT5OjZ0Yk5/NsuBwcHvrOJpOnbrd/PWfg2cXCLCQ94YvlXa7NBb4B/dWEuPKbysrNSHv6eRFfkVys2XWQf/MPkbQdTA6q/eGmczqBk9IU7N37at08fK+UGycyNGzcGjh5fZetuX1OctHtbcHCwRV4Wt0+0a6lxYS18dYFdOHFI/TMbwc3PJmnt6vtdXnx+enteU+BwfXYbWKrQOnSL0cxMBhsHsnl2/xFjFz05jPchegMLP15j33x//SU9oYkvQMn1yJ+eO314b/tfGb9pbSD5Qnv6xTmbbAbTHiPgwqFHC7dv/+pTi7wshz3CvXv3/vTTT56entOnT/fz87s9wZUrV7788suamprHHnssPj7e/GR1dfXGjRtzcnJ8fHymTJkSGBgIAJcuXdqxY0fDHz722GMhISHc5bz1vO3g9Z7MF3aGS3ZOAGC0c62sruY7U0iyTJSC0hYAHJxd53apHdqB/3lqFQOPBjPZQfVLLvsAANg519bV850pJFlVNbXU1RUAwN61qrrGUi/L1Rzh1q1bJ0+e3LVr14KCgr59++r1+iYJNBpN7969TSZTQEDAyJEjDxw4YH5+9OjRu3bt6tOnT2lpaUxMzM2bNwEgJydn/fr1pX8zGAwcZbttlr4+2+vT0fbbZtsf+vTpyRP5zg6SrMeemMSsfdBj24sh2tQHHniA7+z8499TJ3dMWee6Yw6z7qH5c1/hOztIshbMeUG15WWnHa+qt05f+lqrRv5ahXIjKipq06ZN5p/79+//6aefNkmwaNGiRx991PzzypUrR44cSSk1rze5fv16w4t88803lNLdu3f369evNe87c+bMdevWtTKTFRUVrUx5Vzdv3tx/8EjHjfqThaylXlOYLFho8mGpQntkr3H+LxePHTtmMBgs8oIWVF1dffjw4Yk/3HjzhNEiL4jftDaQfKF9d8kU823hocOHS0pKLPiynPQI9Xr96dOnhw4dan44dOjQQ4cONUlz+PDhxgkOHz4MAC4uLkFBQcePHweA/Pz8mzdvRkREmNPodLrly5d//PHHeXl5XOS5nfz9/Qff3++tOIc3TuCqOcSJE4X0eCFd+EBYnz59lErBLXOzt7fv37//ymH+n+WwN3C7LeKAgYVF6ezKQZ4D+vd3d3e34CtzUp3y8/MBwMvLy/zQx8fn4MGDTdIUFBQ0TlBZWVlRUeHi4rJr164RI0a88sorRUVFH3zwQWRkJAA4OzsnJCSwLHvkyJH58+f//PPP999/f7NvbTAYvvvuu9OnT5sfxsbGTp06taV81tTUKBSK9n3WW0wKgDWZyl+u1Cb6SvbQYosXmhxYpNBeO6Z4K4JCfX21gOfg3AlMCWGWnmTXxrX3ihC/aW0g7UL77AIT4kT6utXf00oMOzs7hrlLl4+TQGhjYwMAJpNJpVIBQH19va1t06OZVCpVw63h6+vrCSE2NjY1NTVPPPHEjBkzJk+enJWVNW3atJiYmLi4uEGDBg0aNMicePHixQsXLrw9spoxDBMQEBAXF2d+GBYWdvtbN2g2Y+20tBddfIYMD1Txv4yBG1wUmuS1v9D23ICCWvpMFyL8c/4WxEC3HezLPRTd3NpVCfCb1gYSLrRKA6zMZn8axtja3lvYIq1YWM1JIPTz82MY5saNG2FhYQBw8+ZNf3//Jmk6dOhgXghjTuDp6WlnZ7dv376SkpIFCxYAQFBQ0OjRo7/77ruGqGbWp0+fr7/+uqW3VigU/fv3f/7551uTT4VCYfGrpwmhsCbLuPMamRAi+BarTbgoNMlrZ6GxFBZmGFf0ZmxVIvhSedrDKxFk2Wm6bUi7vif4TWsDCRfa2rPsYH+mlzcnn46TemVraztixIitW7cCQG1t7Y8//jh27FgAqKqq2rNnj7kjOGbMmO+//958DsXWrVvNCby9vcvLy7Varfl1Lly4oFarzX9ofoZSumvXrp49e3KRbYsgACviFAvSWYNkB0eRtW2+zKoYGBsogihoNrsHc0xHj+lwphBZRlEtrM8yLY3hqgpwNeW+dOnSkSNHnjlz5uLFi8HBwaNHjwaAvLy8kSNHlpaWurm5PfXUUxs2bLj//vu9vb2PHz9uXiwTGRn55JNP9u7de+jQodnZ2eXl5f/+978BYNq0adeuXQsICMjNza2srPz11185yrZFDPQjwU7w5Xn2+W6iabmQYNWzsDiD/WKAQkSD7XYKWBjNzD9hOjBacIt6kBgtP216MpQJdeGqEnB4skxxcfHhw4c9PDz69+9vnqusq6vLycmJiIgwd94NBsPBgwdramoGDhzo4uLS8Ifnzp27ePGit7d33759zavjqqqq0tPTtVqtv79/XFyceQ6yWVY+WaYlZ0roiN+MFyaonLk59ZFHkj+6ggvtKbR1Wezem+xPD4gsopgo9NxhXB2veKCtu/7xm9YGkiy0P/U0dpcx61GV2p6rt8Aj1rj60jyZbOruTt6KklqnUJI1jWttLrRKA3Tebvh1uDLKU0Qdwr/s/JNdlsFmPKxk2pR3/Ka1gSQLbcoBU5gLWczZuCjg3Se4804sszrTpLPYGUBIjv6TaRrWgRFjFASAR4IYByVsu4Kz5ajtMkvo3pvsKxHchioMhFwJdiZPhDLvn8X99aiNCmvho2yW0wthrq3orXjzJFuPoRC11fwTpgXRCq7nmERcx4RvcYxi4wX2il5Sg8/IapZlmKaEMSHOouwOmt3vS7q6wRe5GAlRWxwuoNll8O8unMcpDIQc8raDGfcxyzKwFUD37E893XyZnR8p+j1h78cplp8y6YV1Tj4Sh/knTO/EMrbcVwIMhNx6rafij5vsmRLsFKJ789ZJdnYPhQ9ny+SsJsKDDOnArD6Hl4Po3vzwJ1tpgCescjIJBkJuOalgXk/FwpPYCqB7cLaEJmnYl3tIpHou78WsO2fS4sIx1GomCgtOsv/XW9G2Jcf3SiI1TcheuI/JKqUH8rFTiFrr9eOmRTGcLxCwmiBnMimMee8MLhxDrfXVBdbLDoZ3tNIEOQZCztkwsKQXM/+ECSMhao2D+fRCOUzjfoGANS2MVnx7ib1cgZUA3V2tCZZlsCt6W2+CXFKVTbAmhTL1JvjpGg6QorugAPNPmN6LY2ykVTW97OCl+xRLcOEYaoV1WWwfHxLvY7310tKqbULFEFgeq5h3gjViO4DuaMdVtsYIjwVLsGK+GsHs17CnirFTiO6krB4+yDQt62XVKiDB+iZMozoRfwf45hJGQtQiE4VF6ex/+lhpgYCVOangzSjFgpM4U4ju5L3TpocDmXbezPJeYSC0nvfiFIvS2Roj3/lAQvXf86zaHoa29ZRq4XuuK3O+HJI02ClEzdNU0y/PswuirR2YMBBaT29vEudNPs7BTiFqRo0R3jnFrogT/Q76O1AxsKwX8wYuHEMtWJzOTu/KdHS09rUgBkKrei+Oef+MqaSO73wg4VmTxSaoSR8rLhDgxcRQxkhh1594OYiaOl9Of7zGzu3Jw7UgBkKr6uJKxgYy/8GTuNGtSutgzTnT29ZdIMALAvBurOINXDiGbvPGCfb1SIWHLQ9vLf2KJzRLezGf57I3qnBwCP3jndOmR4OZcFeJdwfNhnckHR3h64sYCdE/jhfSE4V0Rjd+QhIGQmvzdyDTujDLT2ErgP5ys4puusgujJby7GATK3orlmSw1bhwDP3tjROmZb0YeyU/746BkAdvRCl+uMbmlGGnEAEALExnn+vG+Ir/fO3Wi/Ui8T7kw2y8HEQAAL9ep5pqmBLGWzzCQMgDNxt4NUKxKB1bAQS5ZfS36+zcCBl1B83eiWX+cxYXjiFgKSw4aXo/jlHyF44wEPJjVncmTUeP6bBTKHfzTrCvRypcbfjOh9V1diWPBDHv40ncsvfdZdZGAQ8F8hmMMBDyw04Bi2KY+SewFZC1NB09W0Jf5GmBAO8Wxyi+PM9ex4VjMlbPwpIMdkWcgt91YjKtgULwTGemsAZ+v4GtgHzNP2FaGmONG3ALk58DTO/KLMOTuGXsk2z2PjcyyI/n9dIYCHmjILA8lnn9uInFUChLP+WxxXUwib8FAkIwL1KxO4/NxoVjslRpgPfPmt6O5b8K8J8DOXs4iHFSwdYreEUsOyyFhSfZ9+IUPA8J8c3VBl7vqVhwEquAHK08a3qgAxPpwX8dwEDIsxVxirdOsnU4Vygzmy6yzjYwuhP/TQDvZtzHpBfRFC12CuWlsBY+ymYXxQgiBgkiE3I2wJd0c4MvzuMVsYzUs/C21M/Xbj07BSzBhWPyszTDNLUzE+IsiGtBDIT8WxGnWH7KpDfwnQ9kLR9msT09SD+1IJoAIXgqnCmpg1+vY6dQLq7q6ZbL7Dw+ztduFgZC/kV4kGEdmFWZ2CmUBb0B/pNpWi6ABQLCoSDwbiwz/wQuHJOLt06ycyIUPoI5TQlroyC83YtZn2XS1vCdD8S998+YRnZkurtjd/AWYwIZFxV8exkvB6XvTAlN1rCzugso+ggoK3IW5EwmhzPvnsZpEonLr4ZPc9jFwlggIDQr4hSL0nHhmPS9nmZaHKNwVvGdj0awQgrFgijFd5fZyxU4NiRly06Znu3CBDhhd7AZ/X1JhDv5NAc7hVJ2MJ9eqoBnuwgr9AgrN3LmZQczuyuW4Ckb0nWxnH5/lZ0XKZQFAgL0f72ZFWdMFbhwTKIowPwTpvfiGBuBRR6BZUfeXunB7Newp4qxUyhNb51kX41QePJxA26x6OpGRnRi/nMWh0el6furrJGFx0IEF3cElyE5c1LBW1GKt05iKyBBJ4toio4KaoGAMC3vxXyawxbgwjHJMVFYnM6+x/f52s3Caiks07syF8vh5wsVpaWlfOcFWYbJZCooKHjjuHFxNOPA0w24RaSDI5kSziw7WVdQUMB3XpDF6HS6T87WdHSEoR0EGAcB66WwqBiITF/78Btb3B1sxg5O+GLtSr5zhNolNzd32PhJlQ6+VaWFXx7YBeDPd45E4EH29NDJ07b5+HvTiqN/7Pbw8OA7R6jtjEbjsHGPZxXoi3Xa5UsWAYznO0fNwB6hsFRXVx/cvdk4/0jhywd+PJZz+fJlvnOE2uXVxe/dGLe+7LldxlFvvvvBWr6zIw5vLF7KPre5+LkfL0ZNW/3x53xnB7XLnj17MkhQ4fO/sPMOfrxqBd/ZaR4GQmExGo1EZQeEAQCwc6qtreU7R6hdauvqwdYJAKiNY1UN/jdbpa6uDmwcAYC1carBQhO5uro6k40jAIDSxmQS6Kp4HBoVFhcXlwcH9N793/El1CHKH+677z6+c4TaZfn8l4dPfro+uK/H9dQFP27hOzvi8N5bc6fOHl/iF+t+M+3lvT/ynR3ULqNGjer4/rqL32g8yi6+8uK/+M5O8wilXC3WLy4uzs3NDQ0N9fX1bTZBfX19RkaGi4tLk+b+ypUrWq02JCRErVY3PGkymTIyMmxtbSMiIghpcbp11qxZ4eHhM2fObE0O9Xq9s7Nz6z6NVZ0/f37Uz9Xbp0bFeAluYlmwhSZYD/5QFFVzbv6YWCcnJ77zIhpFRUVzf8z27hKzsj8W2j0QZvX85kL9V3tPbRrboWPHjnznpXlc9Qh37Njx3HPPRUVFnT59+p133nnuueeaJLh27drgwYN9fX21Wm3Pnj23bdumVCpNJtPEiROPHz8eERGRlpY2c+bMhQsXAkBRUdHgwYNtbGyqq6t9fX1/+eUXe3vBHNfKgS5duowsNiXlUwEGQnRPWAppVa5rhmMUvDdeXl6PD+q1NAs3XUrBIa3ikcS4jh2FOxPHSc4MBsPMmTM3bty4b9++P/74Y+7cueXl5U3SLF++fMiQIUePHj179mx2dvbu3bsBIDk5OTk5+dy5cz/99FNSUtLixYvNuwg++OCD8PDwkydPnjlzprq6etOmTVxkW1AS/UiyRqDj6aj1ThdTHzuitsNDEu5Zby82u5SW1fOdD9RuSRqa6Cfoa3pOAuHRo0dZlh01ahQAxMTEhIaG/vrrr03SbNu27dlnnwUABweHxx9/fPv27QBACLG3tzf39jw8PBQKBcuy5sRPP/00AKhUqsmTJ5sTS1uiP3NUSw0YCkUuKZ8O9hd0EyBYNgz08SFHCvAaQtzyKmmlkd4n7NutcDI0ev369aCgoIaZvODg4Ly8vMYJysvLKyoqgoKCzA+DgoL2798PAIMHDx43btyYMWOio6MPHDiwbt06T09PSunNmzcbJ75+/XpLb82y7Pnz5/ft22d+6Ovr26NHDwt/PKvwsIVgZ3KyiMb7CPoLhO4sWcNOE9j5wiKS6Mck57MPBuDprCK2X0MH+zMCb8U4CYTV1dU2NjYND21tbaurq5skAICGNHZ2dlVVVQCg1+tzc3NdXV0dHByUSuXZs2cppUajsb6+/vbEzTIajUlJSdnZ2eaHAwcObIigt6uqqrrDuhveDfBW/v6nMcLByHdGbiHwQhMUIwtHCmw+jq3DQmuDqqqqvu6KOSeVld3xvLXWEuA3bW+eKsGbrazk7eRIBwcHhrnLxSgngdDX17e4uLjhYXFxsZ+fX+MEPj4+CoWipKTEfGZEUVGROcF///vf+vp688jnK6+80rFjxyeffHLAgAGenp4lJSXmv21I3CwbG5uZM2e2ctUopVTISxgeCKBrs0xLnOz4zsgtBF5ogpKqoyEupkBPR72exUK7V5TS+30crh001CqdvIRVCYRLgNXzaJFxWW8bJ2HfeoyTQZuYmJhLly4VFhYCQF1d3fHjx3v37t04gUKhiImJOXz4sPnhkSNHzAkqKyvd3d3NT9ra2jo6Our1egDo3bv3kSNHmiSWvIF+JE1Ha/EIbtFK1uAEYbsoGeinJocKcKpcrC6WU0ohzEXotYCTHmGnTp3Gjx8/efLkl19+edOmTTExMTExMQDw1Vdfff311wcPHgSAV1555ZVXXnFxcfnzzz//+OOPDz74AADGjh27YsWKDz74IDo6eufOnSaTqV+/fgAwZ86cxx9/vEOHDhUVFd98801qaioX2RYaZxXc507SdHSgsBdcoZYk57Mvd8f5rXZJ9GeSNfSRIL7zgdpELIvFuJrG//LLLwcMGLBhw4aQkJBdu3aZn+zWrdvDDz9s/vmJJ55Yu3bttm3bcnNzDxw4YN5oGRkZefjw4atXr37yySfOzs6pqamurq4AMHTo0P/9738///xzWlranj175HPeSqIfSc7Hy2FRqmchTUf7+4qgFRCyRD+SpMGFo2KVrKGJYgiEHJ4swwtpnCzT4I+bdPkp06EHBXQSnvALTSAO5tN5J0zHxigBC61NzIXGUlB/a8gcr/KV8hEaFiOobxoF8PvWkDZWGSjsCULAQ7cFrr+anCqmVcJaN4paJTmfFfgmYlFgCPRXM3i4hBhllVInFRF+FAQMhALnoIRoT5KilVSvXSaSNDTRH+uXBST6k+R8rALiI6LFYlhRhW6wP561Jj7VRjhVTPupxdEKCNxgf5KM04QilJwv9JPVGmAgFLpEPyYJL4fF5qiWRnsSRwHN7YpYd3eiN9BrlVgLxISlcCifFcuKdwyEQhevJnj0sOgka1ixDAoJHwEY5MccwMtBUTldTNX2xN9BHLUAA6HQ2TDQ2xuPHhaZpHya6IeVy2IScXRUbJLyxbFxwgzrqggk+jO4m1BE9AbILqV98LR0y8HdhKKTrBHTqmkMhCKAiwXE5WA+7e1N7PBIGcvp7EoowKUKrAXiYGQhRUsHimdQRDQZlbM4L3JFT4tq+c4Hap3kfBY3Tlhcoh9eDorGiSIa5ExEdFQ6VlcRwKOHxUVE26dEBHcTiojoqgAGQnEwHz3Mdy7Q3ZXUweUKGuslplZAFIb4k/0aFuuAKCTns+JaLCamvMoZLhYQiwP5bH9fosKKZWkBTsRJSXLKsBYInRiPm8f6Kg7RnkRXSwvwTt2Cl4wnq3Em0R8vB0UgVUvvcyduNnzn415gjRUHPHpYLJI0ojlWSnRwvYwoiPG4eQyEooGLBYRPVwP5NTTKU2StgFgM9mcO5rM4TyhwYjxuXmTZlTPcTSh8+zXsQF9GgXGQG34O4G1PzpRgLRCuaiOcFuFx8xgIRQOPHha+ZFEdKyVGeDkocEe0NEqEx81jIBQNPHpY+ES3fUp0Ev0IHjcoZCI9bh4DoZjg0cNCdr2KVhhod3fxtQIikujPHC6gRgyFQpUszuPmxZdjOcPdhEKWpKGJfgyGQU552kKQE0kvwlogROI9bh4DoZjg0cNClqzBCUJrSPQneKtqYTqYT/v4iPK4eQyEIoOdQsE6kI87CK0h0Y/ghlphEt3Jag1EmWk5w92EwnSpghpY6OyKgZBzA/2YYzpaZ+I7H+g2SaIdFMFAKDJD/EkSHj0sPMkaOkScTYDouNpAVzdyvBArgbCU1MEV0R43j4FQZPDoYWHCHYTWNBgPHRUeUR83L85cyxu2AkJDAZI14jtfUbwS/RjcTSg0oj5uXqz5ljPcTSg0OWXUQUmCnDEQWskAX5JeRKuNfOcDNZKkoYNFey2IgVB8Ev2YA3j0sJCId42ASDkoIdKDpGixDgiF+bj5SNEeN4+BUHz8HMAHjx4WkmS89ZLVDfbHs9YEROzHzWMgFCU8elg4WAqHCthBGAitK9GfwSogHGJfLIaBUJTw6GHhOFNCvexIB0cRtwJilOBDzpXSCgPf+UAAIP7j5jEQipL56GEDhkIBEHsTIFK2CojzJkcKsFPIPwkcN4+BUJTMRw9n4NHDApCcjxsn+JHox+BZa0IggePmMRCKFR49LAQmCke1dKA4z1cUO6wCAiGB4+axAosVHj0sBCcLaYAj8bbjOx+y1MebXK6gxXV850P2JHDcPAZCscKjh4UgSeSL5URNyUC8DzmEq8Z4JY3j5jEQipX56OE0PHqYV3iyGr8S/Rm8GQu/kiRx3DwGQhHD3YT8qmchTUfvxwlC/mAV4J3YdxCaYR0WMTx6mF9pOtrFjbjZ8J0PGYv2JJpqqq3hOx9yRQEOSGJQBAOhiOHRw/xKwh2EfFMQGODLHMDLQZ5I5rh5rgJhXV3dyy+/HBoaGhsb+9NPPzWbZtOmTZGRkeHh4QsXLmRZFgBu3Lgx7Fa//vorABw7dqzxk8ePH+co2+LioIQoTzx6mDfJ+WwijovyLdEPR0d5I5lrQSVHr7t8+fKMjIwDBw5kZmY+8cQTZ86cCQ4Obpzg+PHjs2fP3r17t5+f37hx4/z8/F588UUPD4958+aZE9y4cWPatGmff/45ABQWFubn569Zs8b8q5CQEI6yLTrms9aGdlDwnRHZqTVBehFNUEuhFRC1RH/yUQ72CPmRrKHjg6VQBTi5nqWUfvHFF0uXLu3UqdOoUaNGjRr11VdfNUnz+eefP/XUU/379w8NDX3zzTc/++wzAHBwcBj6txs3bgwePLghfLq5uTX8ysvLi1EA1n8AACAASURBVItsi1GiP4M36eXFkQIa6UGcVXznQ/YiPEhFPb1RhbXA2szHzQ/0xUDYgqKiIq1W26tXL/PDmJiY7OzsJmmysrJiYmIaEuTk5JhHR80opV9//fWzzz7b8Exubm7fvn1HjRq1ceNGSvFL/5cEH5KFRw/zITmflcagkNgRgPt9cRMFD6R03DwnQ6NFRUWEEGdnZ/NDNzc3nU7XJE1xcbGLi4v5Z1dXV4PBUF5e7u7ubn4mOTm5uLh43Lhx5odhYWGffPJJaGhobm7u3LlzS0tLZ8+e3exb19bWzps3b9GiReaH48aNW79+fUv5rKqqIkT0/8UYD9Xeq9XD/a00OiSNQmu/fddViyNNlZWtKnYstDZofaEleCr+uEYe9sVlY1b9pu25qhjgTSora63zdm3m4ODAMHfp8nESCN3d3SmllZWV5lBXUVHh6el5exq9Xm/+Wa/XKxSKhrgIABs2bJg8ebK9vb35Ybdu3bp16wYAMTExdXV1a9asaSkQ2tnZLV68ePr06eaHNjY2jo6OLeWTUurk5NTGDykYQzuyx0qV4ztbaZpQGoXWTpUGyK4wJAbY2reuAmGhtUHrC21kMF2da3JywpPurPpNSykxTg1npFHsnAyNent7u7i45Obmmh/m5OTcvrzF3L1rSBAcHKxQ/NWUl5eX//DDD43HRRtzc3OrqbnTviEHBwf3v90hCkrGYDx62OoOFdDe3qSVURBxrYsrMVG4osdaYD1GFo4U0EFSWTXNycdQKBSTJk1auXKl0Wi8ePHiDz/88NRTTwGAVqt96aWXzGFs6tSpmzZtys/PN/fwzAnMvv3227CwsOjo6IZnDh8+XFRUBADXrl177733RowYwUW2Rao3Hj1sdcka3DghLIP8CK4as6b0IhroRLyk0BsE4G4f4fLlyysrK728vPr27bt48eLIyEgA0Ov1u3btMhgMADB8+PCnn366a9euarXax8dn7ty5DX+bnJw8a9asxq+2d+/ekJAQR0fH6OjomJiYd999l6Nsi5GSgQQ1Hj1sVXjWttDgbkIrS8qXyA5CM8LpCsy6ujpbW9s7JGBZ1mg02ti06pSq2tpaO7u7XIHMmjUrPDx85syZrXlBvV7fsKJH1FaeZa9X0XXx1pgmlEyhtVlZPQRsNhRNUdm0+jISC60N7qnQruppv5+MN59USadtbhOrfdMe+M04szvzUIBExkW4/Rh3joIAwDBMK6MgANw1CspWoj+OC1nPgXy2n5q0PgoiKwh2JrYKcr4Ma4E1mI+bH+ArnTognU8iZ9GeJL+aFuDRw1aRrKGJ/lhxBGcwThNayzHJHTeP9VkKzEcPH8RpQquQzPmKEpPoT3BbvXUkS64KYCCUCFwsYB26GrhRRaM9JdUKSMNgf5KsYVmsBNyT3nHzkvowcoa7Ca0jOZ+9349RYBwUHn8H4mlHMkuxFnCr2gjpRbS/JI4YbYCBUCJ6eJCKeppXia0At5I1VAK3IZUqvGG9FaRoaZQncZTWaRIYCCXCfPTwAewUcixZWtunJCbRD6cJOZecL4Vb0jeBgVA6cLEA1zTVtLiW9nCXWisgGYn+zKF81oSVgEtJUlw1LbXPI2eDcTchx/Zr6GB/hsE4KFTedtDJiWQUYS3git4AWaW0r7fU6gAGQukwHz18uQJbAa4ka/BkNaFLxN2EXDpUQOOkeNw8BkJJwTkSTiXn40oZoUv0J8m4oZYzUj1uXoIfSc4ScdUcZ67oaZ2JdnXDQChog/yYVC2tx1DIDakuFsNAKClD/Ml+DW4p5kSShg6W3BoB6XGzgXBXclyHlcDyyurhcgWNk9wEIWAglJhAJ+KgJLl49DAHcAehWAzG5dPcSNawCWqikmLQkOJnkjdcLMCRgwXSHBSSnkQ/JlmDY6OWl5xPJTlBCBgIpQd3E3Iht4wqCAQ7YyAUgQG+5EQRrTHynQ/JSZLuqmkMhFIzxJ85gEcPW1qShg6VaBMgPU4q6OlBUnGa0KJ0NXBTusfNYyCUGj8H8MKjhy0tOV+y18KSlOiHmygsLDmfHSjd4+YxEEoQHjFjWRTgYD47CFfKiEeiP4NVwLKkvVgMA6EE4W5CyzpbQt1tSUdHybYC0tNPTTJLqN7Adz4kJEmiOwjNMBBK0GB/5nABHj1sMdK7H7fk2Smglxc5qsU6YBmaalpaR7tL97h5DIQS5GkLnZxIOh49bCF4spoYJfrjJgqL2XeTJvpJ+bh5DITShHcotRQThcMF7ECJbp+SsMG4odZyJL9YDKu3NOGqOUvJKKIdHYnanu98oHvU24dcKKeldXznQxIOSHqCEDAQStUgPyZFS+tMfOdD/JIkvVhOwmwY6OtDDhXg5WB7Xa6gdSbaxVXKtQADoTS52kAXV3KiEIeG2is5n5X2oJCEJfozeMpS+yXn0yFSP25e4h9PzhL9SRK2Au1jYCFVSwf6YjURJZwptwg53I8aa7hk4dHD7Zemo+GuxN2W73ygNunlRa5XUW0N3/kQuQMyWDWNgVCyBviSdDx6uH2kehtSmVAQ6K9mcJqwPXLKqIqR/nHzGAgly0kFER4kBY8ebodkDSvV+87IBJ6y1E5JGjpEBteCWMmlbLA/wdHRNqs1wcki2s9X+q2AhOHtOdtJ8jsIzTAQSlmiH66aa7sULY3wIC4qvvOB2iHSk5TU0ZtVWAvaggIcymclP0EIGAilLQGPHm6HZHk0AdJGAO7Hy8G2OltCPWxJBxkcN4+BUMrsFBDrRY4UYCvQFkkamij17VNykOiH04RtlCSb4+axnktcoj+DZ621QZURMktogo8sWgFpG+xP9mMgbBM57CA0w0AocbZZez6dPXHqC7N1Oh3feRGNLdt39B31qPOPb9RVlvGdF9ReDuV5uv/O6PfQE0eOHuU7L6Jx4cKFsZOm7Vn6lHvBab7zYg13CYRVVVUXLlwoLy+3Tm6QZeXk5Pzf/63QP7j0f0z/MZOm8Z0dcUhLS5ux4vNzie8WePR44l8z+M4Oaq8Hxk+qiXw4Je6N8dPnFBQU8J0dEWBZdujDE3d3esowdM6UZ/+l1+v5zhHnWgyEX3/9defOnZ2cnLp06eLm5hYQELBmzRpKcYRBTDIyMvQ9xoJPGBs99trNfL6zIw6paSdKIx8H72C276RzORf4zg5ql7q6urI6Ct2GgP999eGJmZmZfOdIBLRabb2zH4QlQKdIY2Dc+fPn+c4R55TNPrty5crXX3+9T58+K1euVKvVxcXFv/zyy5w5c/Ly8latWmXlLKI2i4uLc3nv6aKw+4kmKyywE9/ZEYcB/eLdN7xWEhitunQkumd3vrOD2sXW1tbLQVl45ifqora9sD8y8lW+cyQCvr6+dlVayNoLtg6qaye7dv2A7xxxj96mrq7OxcXlueeea/L822+/rVAodDrd7X8iHDNnzly3bl0rE1dUVHCaGSFISj7QZ9xU/8feKCkpscgLyqHQVvzvV7eBU+YuWKrX6y3ygnIoNIuzVKFpNJpnXnpFNfDZP1LSLfKCQmapQrt69ardsBmjJ0/PysqyyAsKXDM9wsLCwoqKihkzms6OvPjiiwsXLrx27Zq3t3dr4uv27dszMjK6des2adIkpbKZN0pPT//xxx+dnJymTp2qVqsBIDc399ChQ43TTJgwwc3NDQCys7O3b9+uUqkmT54cEBDQ+kgvc4mDBu7tN9D/O4OjK+4Mby3aY/gzKx5Y2VfBd0aQBfj5+W1Y/0HB78Yq3AzTakaPQO9n1vw8sfkhQ+lp5pvh5ORECCkuLm7yfFFRkUKhCAwMbM3rzp07d/ny5Wq1+rPPPps6dertCZKSkoYMGWJra3vp0qXY2NjS0lLzW6T/bfv27a+++qo5gmZkZMTHx7MsW1hY2KtXrxs3btzzB5UxZxWEOpPTxTi/21qpOhqvlsWqcfmI92FStVgFWitFSxNkVQWa7ScOHTq0R48ely5danjm5s2bAwcOnDVrVmu6mUVFRfb29uY/Lykpsbe3v3jx4u1vsWrVKvPPDzzwwAcffNAkwfTp05955hnzzxMnTnzjjTfMP0+aNKnh59vh0GizXjhiXJNpsshLSb7QWEq9v6m/Xsla8DUlX2hcsGyh7bvJ9v/JYMEXFCZLFdrzlmsxRKH5sYKvvvqqvr6+c+fOMTExo0aN6t27d3BwcHp6ukajmfC3pKSkloLrsWPHOnbsGBoaCgDu7u5xcXHJycmNE7Ase/DgweHDh5sfDh8+/MCBA40T1NTUbNu2bdq0v1b8JycnjxgxoiFxk1dDdxWvJql4G4rWuVRO7ZWkowyOlZKVvj7kdDGtx7MlWidVZj3C5oeAO3bsmJ6e/s033xw4cODGjRssy0ZERADA1atXG9JUVFS09KL5+fk+Pj4ND9VqdX7+LWv3i4qKDAZDQxq1Wq3RaBon2LZtm4+PT0JCAgAYjUadTtcwMXl74sYMBsN333135swZ88PY2Ninn366pcS1tbUqlSxmzmJcYX6BorbWAqeOSr7QDt4kfbxIbW2tBV9T8oXGBcsWmgIgxIlJ09TFeUn5itAihaY3wGW9ootjnUUrAW9sbGwY5i7Twy3OhTo5Ob3wwgsvvPBCG95YqVSy7D+XXiaTqcliGfPDhjRGo7HJP2/Dhg3Tpk0jhAAAwzAKheIOiRtjGCYgIKBXr17mhyEhIQpFi0seFArFHX4rJZ3dgKVQUKfo4NDel5J8oR0vgngfsOxnlHyhccHihZaghuPFpK/agi8pOBYptHQtRHuAvUoi31hzHLkzThYF+fn53bx5s+GhRqN56KGHGidwd3e3s7PTaDTmTqFGo/Hz82v47cWLF1NTU7ds2WJ+yDCMuRfYvXv32xM3oVAo+vfv38r4rVKp5HOd3sfHdLyYTHBt78I5yRdaWpFx+n0KlcqS40KSLzQuWLzQ+vmyP+VRlVTa92ZZpNCOl7AJvhIvqCY4WU/cv3//8vLykydPAkBeXt7p06fN04EFBQXmQUtCyOjRo3fs2AEALMvu2rXrwQcfbPjzDRs2jBw5snG0e/DBB3fu3Gn+eefOnY0To1aK98FpwrurMMAVPe3pIaPZEfmI98E7sbRKqpaNl9lx85z0CB0dHZcsWTJ27NixY8f+/vvvs2fPNke1H3744bPPPjt9+jQALFiwYMiQITdu3MjLyzMajRMnTjT/rdFo3LRp0yeffNL4BV977bWEhITy8nK9Xn/58uX//e9/XGRb2hLUZG4aLhW4izQd7eVFbHC/mRSFuhCW0utVtBOuhGoZBUgrpF8PlFcd4Gq/5OzZswcOHHjq1KkpU6bEx8ebnxw3blyfPn3MP0dFRWVlZe3du9fZ2XnkyJG2trbm5+vq6jZu3Dho0KDGrxYSEpKdnb1nzx4bG5sRI0Y4OTlxlG0Ji/UiWaW02ggOctkj2xYpWiq3a2FZ6evDpGjp4yH4L25Rdin1tCVqe77zYV0cNorR0dHR0dGNn/Hz82s84Onr6ztlypQmf+Xo6Dh06NDbX83Dw+PJJ5/kIp8yYa+E7u4kvYgO8MVWoEWpOvbFbvK6FpYV8z6ix0P4zoeAyW4rPQDg/QhlJQF3E94RBTheSOPVWCkkK8GHpOD5MneUqpPjoAjWeRmJ9yF4ytQdZJVSLzvibcd3PhBnYr1JdimtNvKdDwFL1WGPEElaPzVJ0eF6mRalammC/K6FZcVOAT08SHoRXg42r7QO8qtpd3fZ1QIMhDLSwZHYMuRyBbYCzcOztuUAR0fvIEVHe3sThfwqAQZCeYlXkxScJmxBCvYIZQDP3b0DGe4gNMNAKC84TdiS4jooqKH3yW9QSG76qclRLYt1oFkpWpkuFpPjZ5YzXDjaklStTAeF5MbfgdgrcIKgGSYK6UW0j7cc6wAGQnmJ9iSXKmiFBe5CITWpOlaGi+XkKUGN04TNOFNMOzkRd1u+88EHDITyomIg2pOcKMRWoKkULY33weogC3jubrPkuXHCDGu+7MTjqrnbGFlIL6K9ZTkoJEMJapwpb4Y8t9KbYSCUnXg1ScXdhLc6U0ID5TooJENRnuSynpbX850PgZHzQbsYCGUnwYc5pqO4bK4xOQ8KyZCKgRicILiVtgYq6mkXN5nWAgyEsuNjDx62JLccW4F/pMr4WliecENtE0e1bF+fVtzKXaIwEMoRHq7RRAqeKSMz8T4kVYsTBP9IlesOQjP5fnI5i8fFAo3kV4O+nnZ2xUAoI/3UTFohThD8I0XeswMYCOUIt9U3lqJl49XyHRSSJy878LQlOWVYCwAA6kxwppjGecm3EmAglKMe7kRTTYtq+c6HMKTqcAehHOG2+gYZxbSLG3FS8Z0P/mD9lyMFgThvchxXzQGAXG/JjXBbfQO8ARkGQpmK98HdhAAAdSY4W0JjZTwoJFs4QdAAb0CGgVCm4n0YHBcCgPQi2lXeg0Ky1cOd5OMEAQDI+0wZMwyEMhWvJicKqVH2fULcSi9bDIHe3iRN9hME1yqpkaXBzrKuBRgIZcrNBjo5ksxSubcCeC0sZ7ibEABStLSfjHcQmsn988sZrpoDgGMYCGUsXo0TBDhBCICBUM7iZb9Y4E89pRSC5D0oJGd9fcjJIrlPEKTIfskoYCCUswQfuZ8vgyeryZybDQQ4kbMl8q0FNUY4X05jZL9qGgOhfHVxIxUGml/Ndz74g2dtowR5n76dVkgj3Imdgu988A0DoXwRgD7e5JiMdxPK/HxFBH+tl5FvIMRV02YYCGUtXs3Idpqwygjny2iMJ7YCsibzbfU4KGKGgVDW4mV8P6YThTTSk9jKflBI5jq7Er2BaqrlWAsowDEd2xcDIQZCmevjQ04X0zoT3/ngAy6WQwBAAPr6kGOy7BReLKeOKtLBEWsBBkJ5c1RCZ1dyqliOrUCqjsUlowgA4n0YeU4TpuC46N8wEMqdPOdIKMAxHe3jja0Agni5LhzFY5UaYCCUO3mumrtQTp1wUAgBAEBvb3JGlhMEuGS0AQZCuYtXk6PyC4Q4QYgaOCqhixvJkNkEQYUB/tTTnh5YCwAwEKIQZ0KB5lXKqxXA8xVRYzI8ZemYjvbyIiqMAACAgRABQF/53ZsQb8mNGpPhubspWhYnCBtgIETmu9XLqBWoMMC1ShwUQv+Q4YbaVC0OivwDAyGS3f2YUrU01oso8buP/hbsTAiBa7KZIGApHC+kfX2wDvwFCwJBrBfJLafVRr7zYS24gxDdro+3jC4Hs0qpjz3xtuM7H4LBYSA0mUxZWVn5+fl3SJOXl5ebm0tp0+9ffn5+ZmZmZWWl+aHBYChtxGAwcJVpWbJVQA93cqJQLq1AipbG47UwupWspgnxuPkmuGoOrl692q1bt0mTJkVFRb300ku3JzAajY8//nhCQsIjjzzSq1evoqIi8/OVlZXjxo3r3r375MmTAwICzHF0z549Pj4+oX87cOAAR9mWLflsq/97UAhbAXSLBDlNE+JZ201wFQgXLFgwfPjw06dPZ2dn79y589ChQ00S7NixIzMz8/z589nZ2WFhYStWrDA/P3PmTADQaDRnzpy5ceOGp6en+fk+ffqU/G3YsGEcZVu25LNe5lwp9bUnXjgohG7Vy4vkltFKeQw24Vb6JjgJhAaDYceOHdOnTwcAT0/P8ePHb926tUmaLVu2TJo0ydHREQCmT5++ZcsWAKioqPj222//85//GI3Guro6BwcHGxubhj8pLCysr6/nIsMoQU1StKwcImGKFpsA1AxbBfT0ICeLpF8JimpBW0O7uWEt+AcngbCgoKCuri40NNT8MCQk5Nq1a03SXLt2LSQkpCGBRqOpr6+/fPmyra3tu+++GxMT4+/vP336dJPpr4OPTp48GR0d7erq+uijj5aWlrb01izLXr9+Pf1vV69e5eDzSZC/A3FUkkvl0m8F8HxF1BKZTBCk6ti+PkSBlaARJRcvWlVVBQC2trbmhw4ODnq9/vY0dnZ/jU/Z29tTSqurq0tKSiorKwMCAi5cuFBaWtq7d++NGzc+++yz/fr10+l0Li4uxcXFjzzyyLx58z7//PNm39poNG7evHnv3r3mh8OGDVu4cGFL+WxYjIMAIM7TJimv3jfoLkcuir3QjhbYzgit1+ut2t6JvdB4Yf1Ci3JmvvtT+VKIiIedWlNoB68rY1xBr6+xQn6EwMHBQaG4y31HOQmEarUaAEpLS728vACguLjY19f39jQNHbuSkhI7OztXV1dzsqeeegoA3N3dx44de+TIkWeffdbDw8Oc0tPT89VXX50zZ05Lb21jY/P666+bJxpbw9nZ+d4+m3QN8GdPlSqfc777nWrFW2hFtVBSb4jt4MRY/XJYvIXGIysX2uAgOjvd6OTsLOrO0l0L7WSZ8c0ohcg/pYVxMjTq7u4eEhKSmppqfpiamtqrV68maWJiYhoSpKSkxMTEEEJCQ0Pd3d0bB0hXV9cmf5ifn+/m5sZFtmUuQQY3o0nRsn28ifWjIBIFfwfibEMuSnqCwMjCqSK8AVlTnPQIAWDWrFmvvfaavb39uXPnjh49+uWXXwLApUuXhgwZkpmZ6eLi8sILL8TFxSUkJKjV6qVLl65atQoA7OzsZsyYMXv27Pfee+/y5cvff/+9ebnp2rVrXV1dAwMDc3JyFi1a9N5773GUbTmL8iDX9LS8Hlxt7p5YpFJ1NF6NOwhRi8xnrXV2lWycOF1CA52JhOt423AYCFUq1apVq9zd3ZOSknx8fADAyclp5MiRKpUKALp06fLzzz9/+OGHNTU1K1asmDBhgvkPlyxZsm7dunfffdfLy+uPP/6IiooCgKCgoO3bt2u1Wj8/v40bN44ePZqjbMuZkoFoL5JWSB/oINlWIEVLF0RjIEQtMu8jeroz3/ngDN6ArFnk9lNdRG3WrFnh4eGtnCPU6/U4c9PYmydMtgqyOOZOoUK8hWZgweMbw42JKutfDou30HjES6GlF9GnD5oyx3PVQ+DaXQttYrJpREcyNRwvB2+BxYH+Ea8mqTqW71xw5XQxDcFBIXRHkR7kWiUtreM7H5zBG5A1CwMh+ke8D3NMR6W6rx5P00B3pWSglxc5IdFt9ZpqWmWkYdKdAW0zDIToH152oLYn2WXSbAVwKz1qDfMpS3znghMpWpqgxkXTzcBAiG4h4TuU4uFqqDXifZhUiVYBvBZsCQZCdAup3oxGU02rjTTUBVsBdBfxapJWSE0SrATmG5BhFWgGBkJ0iwQfIsnL4aM4KIRax9MW1PYku1RqtaDOBJklNBa30jcHAyG6RXd3oq2hhbV858PS8AZsqPUkecrSySLazY04inVjCLcwEKJbMAR6e5M0ybUCeEtu1HrxUhwXwVXTd4CBEDUlvd2EdSbIKqWxXtgKoFaR5P2YcFDkDjAQoqbifRiJLRw9UUi7uREHHBRCrXOfGymsldoEQaqOjcceYQswEKKm4tUkvYgaJdQnxEEhdE/MEwTHJDQuclVPCZBAJ6wFzcNAiJpyUUGgEzlTIp1OIW6fQvdKYrsJcRPtnWEgRM1IUEtqW32qlsVAiO5JvLQWjuK14J1hIETNMN+Mhu9cWMYVPWUICcBBIXQv+vqQ9CJqkMrgKPYI7wwDIWqGlFbNpWhpP2wC0D1yUUGws0QmCKqMcKGcRntiLWgRBkLUjHBXUmWgN6uk0Aqk6igulkNtkCCVc3fTdDTKk9gq+M6HgGEgRM0gAH19mGOS6BTiLblR20jm3F1cNX1XGAhR86TRClQZ4VIFjcat9OjeSebcXVwsdlcYCFHzpHE/pmM6GuVJbPBrju5dmCupNtIbIp8goABphbSvD9aBO8HSQc3r400yS2mtie98tE8qjouitiIA8WrRTxCcL6MuKuLnwHc+hA0DIWqevRK6uJKMInG3AnisFGoPCewjSsHFYq2AgRC1SOw3o6EAx3Q4KITaTgITBHjWdmtgG4FaJPab0eSWUXdb4mvPdz6QaPX2JudEPkGAW+lbAwMhalGCmhzVivhojRS8FkbtY6+Erq4kXbQTBGX1cL2KRrhjLbgLDISoRYFORMmQq3qxtgJ4viJqP1Gfu5uqpXHeRInN/N1gCaE7EfViARwUQu0n6iqQqmOxCrQGBkJ0J+KdJiyrhxtVtAcOCqH2SVCTFNFOEKTqaDwuFmsFLCN0J+K9GU2qlvbGQSHUbgFORCXOCQKWwolC2tsbrwXvDtsJdCcxnuRCOa008J2Pe4eDQshSRLqJIrOU+jsQLzu+8yEGGAjRndgqoKcHOSHCVXMpWhwUQpYh0nN3cdV062FLge5CjEcPmyicLKJ9sBVAliDSHmGqFs+UaS0MhOgu4tUkVSeyxQKZJdTfgXjY8p0PJAm9vMjFCvFNEODdl1oPAyG6i35qJlVLxXU9jBsnkAWpGIj0IMcLxVQJdDVQXEe7umItaBUMhOgu1PbgYkMulIupFcCt9MiyEsQ2TZiqY/v6EAYrQetgIER3J7o5EhwUQpYV70NSRbWbEHcQ3hMsKXR34lo1p6uBkjraBQeFkOUkqJkUHWVFUwlwyei9wUCI7i5BVD3CFB0bj4NCyKLU9uAungkCAwunimlvDISthoEQ3V2kJ8mrpKV1fOejdVJxByHigIhOWTpVTEOdiYuK73yIB7YX6O4UBGK9SJpIVs3hLbkRF0R07i6umr5XGAhRqySoxbFYwMDC6WIah+crIksT0f2YUvFa8B5xFQi1Wu2UKVMiIiLGjx9/9erV2xOwLLt8+fKoqKj777//559/bni+vr7+nXfeiYuLi46Ofv311xueX7t2ba9eveLj4zdv3sxRntEdxKsZUayXySiiYS44KIQsr6cHuVFFS8QwQZCqpQk4QXgvlBy97uTJkwMDA3fu3Llhw4YxY8acPXuWkFv+MR999NHmzZs3b978559/Tpo0KS0trWvXrgDwr3/9Ky8vb9WqVW5ubufOnTMn3rJly+rVq7///vvy8vIJEyYEBQXFx8dzj8B0cgAAHFpJREFUlHPUrL4+5EkdNQk+FKbgxgnEDQWBWG+SpqMjOwn6C3azitaxNNRF0JkUGk56hOfPnz98+PDq1avDw8OXL1+u0+kOHTrUJM3HH3+8ZMmSnj17jhkz5tFHH/3iiy8A4OzZsz/88MOOHTsGDBgQERExceJEc+KPPvpo3rx5sbGxQ4YM+fe///3JJ59wkW10B5624O9IskqFHglTcdU44kyCGI4bPKqlCbhY7B5xUl5ZWVldunRxdnYGAIVCERMTk5mZ2ThBfX39+fPn4+LizA/j4uLMCY4fPx4XF/fdd9+NHj363//+95UrV8wJzp07Fxsba/45Nja2yash6xDFtnrcSo+4E+8jggkCnCBsA06GRnU6nZubW8NDd3d3rVbbOEFRURGl1NXV1fzQzc1Np9MBwPXr148dO9a9e/elS5fu2rWrf//+ubm5dnZ2ZWVlDYnd3d3NiZtVV1e3dOnSVatWmR+OHj16xYoVLSWuqqpqMmCL7iDGVXFYw4xzF26hXa+COpOND6mtrOQ7K7fCb1obCLDQejqS4zpVub5SIax8/aOqquqwxubdaFNlpdB7rlbj4ODAMHfp8nESCF1dXauqqhoe6vV6d3f3xgnMYbKqqsr8fGVlpfkHFxcXOzu71atXK5XK2NjYLVu27N+//+GHH3ZwcGh4wYbEzbKxsXnppZeeeuop80MnJycnJ6eWElNK7/Bb1ERiAF173uTo6CjYQjurY/v7CvF/it+0NhBgoTkBdHA0/mlwjPQQaCSsNtDcCqZ/J1sHrpZ/SBMnQ6PBwcFXrlwxGo3mhxcuXAgODm6cwMHBwcfH58KFCw0JgoKCACAkJMTZ2Vmp/Ot/6ObmVllZCQBBQUEXL15skrhZhBBPT8+Qv/n4+Fj2o8nZfe6kuI7qagXaBACetY24J/AJglOlTHd3glHwXnESCPv06ePl5bVx40YA+O2334qLi0eOHAkAaWlp69evN6eZPHnyunXrWJbV6XTffffd5MmTAWDUqFG1tbX79u0DgPT09JycnISEBHPijz76yGAwVFRUbNiwwZwYWRkB6ONNThQLN9Lg+YqIa/FqQW+rTyticI68DTgJhISQjRs3Llu2LCAg4Jlnnvnmm2/s7OwAICsra9u2beY0CxcurKqq8vX17dKly+TJk4cMGQIAtra233777bPPPhsWFvbggw9+/vnnoaGhADB79mwPDw9/f//AwMCBAwdOmDCBi2yjuwrRZ33x5YbTp0/znZGmampqNmz69uyeLd2dxLDPC4lWD9uy37d+tWvXLpYV3CRccnLy999+HWbI4zsj4kMoZ7dcZVm2pKTEzc2tYajzdqWlpXZ2dvb29o2fNJlMZWVlnp6eTRKXl5crlUpHR8c7vOmsWbPCw8NnzpzZmhzq9Xrz0lbUGr/v3Tdh9pKK2CmeZ7Z8seyVh8c8xHeO/sKybGS/wRc6DDGa2J6FhzIO7xPaIgv8prWBAAutoqIiImHw9YgnHcvzEl3Ldm/+mu8c/WPxuyvX/ZZeFnS/Z9oXx37dHhYWxneOxITDsWSGYby8vO6cptllLwqF4vYoCAANC0cRLz7euKXi0TXQKbK4y8B1Xy4RTiC8fPmyzsa3fsR8AMjfmHnt2rU7zCIj1GYpKSll4Q/QITMrAdJW9zcajXe4yreyTVt3lM1IBkZZ4uC6ZceuBfPm8p0jMcF9l6i1wgI7qPLSAYC5djI0sCPf2fmHt7c3q7sMhhqoq2ILrzR7FYVQ+3Xo0EGlOQusESq0NlRAURAA1D7ecOMcADjcOBkSIKDqKQoC+kcigVsyf+7pydPOrvys3Dlgya9f8Z2df7i5uQ175tVdHyR62TErlr0ltPE0JBkRERGzHhn86ep+xSbbme+u4Ts7t3j/g1VDJz7vYSgeMWTgE4/jKop7w+EcIS9wjpBrer1++kmHGE/yWk+hDCfUmiBsm/GX4QrB7u7Cb1obCLnQdlxl3zvDnhinFM4XbtohU4ATeSW8SrCFJmRCacuQiCyKZv6TadIb+M7H3z7OZhN8iGCjIJKeR4IZIwu/XRdKL+JaJf0pj53VHdvzNsKCQ/esmxsZ4s98kiOI5eO1Jlh1jl0Qjd9kZD0EYFEMsyhdKLdjWZrBzrhP4W7Ldz5EC5sP1BZLYphVwugUfpTN9lOTntgdRNb1cBBjZOFXAXQKL1dgd7C9sOxQW3R2JUP8mY+zee4U1hhh9Tn2rSj8GiNrM3cKFwugU/jOaXZmd+wOtgu2IKiNlsQwq8/x3Cn8KAe7g4g3DwcxJgq/5PEZCrE7aBFYfKiNwl3J0A7MR/x1CquN8MFZ0wLsDiKeEIBF0cziDD47hctPs7O6K9xs+MuBJGAjgtpucTSzhr9O4cc57EA/JgK7g4g/44IYBnjrFF6uoD/nsTOxO9huWIKo7cydwg/56BRWGeGDsyacHUT8IgBv8dcpfPsU+zJ2By0B2xHULoujmbV8dAo/zmYH+WN3EPFvbCDDAPycZ+3LwcsV9Jfr7EvYHbQELETULuGuZJjVO4VVRliVaXozEr+9iH8EYEE0szidtXKncNkpdnYP7A5aBjYlqL0Wx1h7pvCjbDYRu4NIMMYGMkoGfrpmvcvBSxX01+vsjPuwAbcMLEfUXmEuZHgHZn2WlVqBKiOszjS9ibODSEgWRDFLMqzXKVyWgd1BS8LWBFnAohhmbZapwiqdwg+z2MH+TA937A4iARkTyCgZ2G2VTuGlCvr7TVwsaklYlMgCwlzIiI7W6BRWGWHNOROeLIoEaEEUs9QqncJlGezL3RUuKu7fSTawQUGWsTCaWXPOVFbP7busz2IH+zPd3LA7iATH3Cn8keNO4cVy+gd2By0NSxNZRpgLGd2J+ZDLTmGVEdaeMy3E7iASqoXRnHcKl51iX+6hcMbuoEVhm4IsZmE0szaLw07h+ix2SAemK3YHkVA9FMDYcNkpvFhO995kX8LFopaGBYosJpTLTqG5O4gniyKBWxitWHiS5ahXuPQUOxu7gxzAZgVZEnedwnVZ7FDsDiLBezCAOCg56RReLKf7buLeQU5gmSJLCnUhDwZYfvlopQHW4WJRJBKLYhSL0i3fKVySgd1BrmDLgixsQRSzztKdwnVZ7LAOTBdX7A4iERjdiTgoYZdFO4UXy+l+DXYHuYLFiiws1IU8FMCss1ynsNIA67NMb2F3EInHYkt3CpdksHOwO8gZbFyQ5b0VxazPMpXWWebV1maxD3TE7iASk1GdiJMKfvjTMpeDF8rpfg37InYHOYMliywv1IWMsVCnsNIAH2bhyaJIfBZHKxZnWKZTuCSDfSUCu4McwvYFcWJRDPNRtgU6hWvOscOxO4hEaGQn4qSCne3uFOaU0f0a9oVu2FZzCAsXcSLQiTwUwKzNMrXnRSoMsC7L9AZ2B5E4LY5WLGl3p3DZKXYudgc5hk0M4sqiGObjbLY9ncK159hRnbA7iMRqZCfibgs72tEpzCmjSdgd5B6WL+JKoBMZE9j2TmGFAdbj7CASubeiFEvb0SlcmsG+1lPhhN1BjmErgzi0MLrtncK159jRAUxn7A4iMRvRse2dwuwympzPPt8VW2nOYREjDgU6kbGBzJpz99wpNHcH34jE7ycSvQVRiiVt2lOI3UGrwYYGcWtBNPNJDltyj53CNefYB7E7iCRheEfiYQffX723TmF2GT2I3UFrwVJG3Ap0IuOCmLX30iksr4cPs0zzsTuIpGJh9D3PFGJ30JqwrUGcWxB1b53CNefYh7A7iCTkgQ7E0w62t7pTmFVKDxfgYlHrwYJGnAtwIg8HtXamsLwePsrG7iCSmoXRimWt7hQuzWDnRigclBznCf0NmxtkDW+1ulO4+pxpTCATjt1BJC3DOhBPO9h25e6dwqxSekTLPo/dQSvCskbWEOBEHgliVt+tU1heD5/ksLh3EEnSwmjFslN37xQuyWBf64ndQavisLDPnj27Z88eNze3iRMnOjs7356gqKho69atNTU1Y8eODQ8PNz/5+++/V1RUmH/29vYeNGgQAGg0mqNHjzb8Yf/+/f38/LjLOeLCW1FM9A/Gl7srvOxaTLP6nOmhACbEGbuDSIKGdSDedrD1CjsxtMVLvaxSmqKlGwcqrJkxxNWl9++//z5o0KCKiorffvstISGhtra2SYKSkpJevXqlpaXpdLq4uLj09HTz83PmzPnkk0+2b9++ffv2AwcOmJ9MT09/8cUXt/+toKCAo2wj7gQ4kceC7zRTiN1BJHkLoxVLMlhTy53CxRns3J4MdgetjKvyfvvtt999993nn3+eZdm4uLjt27dPmTKlcYIvvviie/fumzZtAgAHB4cVK1Zs37694W/79evX5AW7dOmybds2jnKLrONNc6ewh8K7uU7hqkzTGOwOIkkb2oH42sO2FjqFWaU0VUs3YXfQ6ji5+q6urj569Ojo0aMBgGGYUaNG7d27t0maffv2jRo1yvzz6NGj//jjj4ZfJSUlbdy48cyZM43Tl5WVbdiwYdeuXeXl5VzkGVlBgBOZENL8nkLsDiKZWBSjWNxCp3BxBvsadgf5wEmRazQaAFCr1eaHvr6+qampt6fx9fU1/+zn51dRUVFZWenk5BQUFHTp0qWLFy/OmTPn6aefXrVqFQCoVCpfX9+0tLTc3Nznn39+z549UVFRzb610WjctWtXXl6e+WHPnj0nTJjQUj7r6upsbGza91llp52FNvc+6P0TeaGzycv2lpbg/86QhzqBv019nYXuay8o+E1rA6kWWn9P8LUj/ztf90TwLc9nlUFKAfNFX7Y9VUCqhdYeNjY2hNxlnImTQGh+V0r/aulYllUomnb2CSGNEwAAwzAA8Ouvv5qfvHz5co8ePaZOnRoZGTlixIgRI0aYn589e/a8efN+//33lt7dwcHB3d3d/LOzs7P5ZZvFMMwdfoua1c5CC3CCx4Lhw1xYFv3PV7OsHr64AEdHA8NIc1wUv2ltIOFCWxAJM46RCcGgbPT53j4Dc3uAk027PrKEC41TnARCf39/QkhBQUFgYCAAFBQUNHT+GqdpWPNSUFDg5ubm4ODQOEFoaGhoaGhubm5kZGTj54cPH/7DDz+09NZKpfKBBx6YOXNma/KpUqlUKjzC6N60v9DeiqFRO42v9FQ1zBSuO2t6OAjC3SU7NYLftDaQcKE9EAAdMo07rzOTwv4KWudK6fEi03eDlar2NckSLjROcXLtYG9vf//99+/evRsATCbTL7/8Yu7P1dXVZWZmmvt/I0aMMCcAgN27d5sTmH9llpeXd+XKlbCwsCbPJycnN+y1QGLUyZE8HsKszvxrprCsHj7LYfEoGSQri2MUS0+xxr8btkXp7LxIxh5nB3nCVcEvWrRo/PjxV65cycnJUSqVjzzyCABcuXLl/9u7+6Co6jUO4D/cBUSZFXch9oWXZbyEaaYIQWnR2mUVEGUQSBpuJCCgoZRp03jt3qHSmWDQoik16o69aDAZKkHyqgSogW9FBjmiyGC7A2y8uMjCLrvn3D/OzF4H8O2O7jns+X7+Oss8Ot95ZjkPe36/s+epp54aGBhwc3NLTU3dt29fXFycTCb79ttv6+rqCCEXL15MT09fsmQJRVElJSXJyclBQUGEkLS0tKGhIR8fnz/++OPs2bOVlZWPKDbYxvZF0xYdMW9ZIPCYTnZfsqxRTvPDZlHgk2UyB/kMUtxB/eNv034foJt76UMqu70iwn3/W6h76Nrb26uqqtzd3WNiYlxcXAght27dqq+vX7FihVAoJITo9fqjR48aDIbo6Ghvb29CiMlkOnXq1OXLl4VC4eLFi4ODg5n/6saNG42NjTqdTiaTqdVq6xLgRNnZ2f7+/vd5aXRoaGjSO/3hLh5W07LOWGbSxvU++ufqJc0xQvsehHin/R/svmmnuunUBnPtEl3275K/eztunv8QLorYfdMekUc4CFmBQfioPaym/edwecbWHc5imchZ0NV41L63uuGd9n+w+6bpdDrf51Y5uLob+zRnjn4ZErjw3v/mXuy+aY8IFmaAHTt37qTeOjmyqUyvCDpy591PAPYqt2DvqGqzYeNRS+qXb/1rJ9txeA2DENhhsViI0IkQQjnNGB21x5sHAe5q1GgiTjMJIcR55qhd3j87dWCXErDjn1uy/l0QSUvnuvdfjo/D7ifgnW1Z6cci1xivVk3rPJe7fw/bcXgNgxDYsSFt3erI5VqtduHChbjzCXhIqVReuXDq0qVL/v67xGIx23F4DYMQWCOXy+VyOdspAFgzY8aM0NBQtlMA1ggBAIDfeD0Ig4KCent72U4xxfj7+xuxsP8gDAbD3Llz2U4xxWi1WnxUelCtra0rVqxgO8WUxOtBqNfrb//yNrgfg4ODbEeYYmiaRtMeFEVRer2e7RRTjMViGRoaYjvFlMTrQQgAAIBBCAAAvGZvu0Y9PDw+//xz63Mt7m5sbCwxMRF79x/I9OnTo6Ki8Myz+0dRlKOjo1qtZjvIVGIymUZHR9G0B2IwGHp7e9G0cfbv3z9nzpy719jbd42OjY3V19eznQIAADghNDT0nt+/am+DEAAA4IHgAhcAAPAaBiEAAPAaBiEAAPAaBiEAAPCaICcnh+0MttDf33/s2LErV674+vpO+jB0s9lcXV19+vRpsVg8a9Ys2yfkoOHh4bKyspaWFoVC4eLiMrHg6tWrNTU1ly9fnj17Np6Lzejs7CwtLe3p6fHz87vLTSZtbW1tbW1KpdKG0biro6OjtLRUp9P5+fk5ODhMWnPx4sXKykqNRuPp6Tl9+nQbJ+Sg9vb2srKy/v5+pVI5adM0Gk1FRUVbW5tYLMav5z3QPNDR0eHp6bl27dqoqKiAgID+/v5xBRaLZfny5U8//XRaWppYLD558iQrOTllYGAgICAgIiIiMTFRKpVeu3ZtXMGnn34qk8kSEhLWrFkjEonKyspYyckp1dXVYrE4LS0tODg4MjKSoqhJy7RarYeHh0QisXE8bvrxxx8lEsn69esDAwNjYmImFlAUlZGR4e3t/corr6xatSo/P9/2IbmmpKREIpGkp6cvWLDg5ZdfnlhQWVnp5uaWlpaWkpLi5uaGc9rd8WIQbty4MSMjgzmOiIjIzc0dV1BRUaFUKg0GA03Te/fuXbp0qa0jck9eXp5arWZO5Rs3bszMzBxX0NXVZTQameP8/PzAwEBbR+Se0NDQ/fv30zRtMBh8fHxqamomLYuNjd22bRsGIWPRokUHDhygaXpoaEgmkzU0NIwrOHjw4Jw5cwYGBlgIx0kURQUEBBQXF9M0PTg4KJFIzp8/P64mOjo6JyeHOd6+fXtcXJytU04pvFgjLCsri4+PZ47j4uLKy8vHFZSXl69cuZK5+hcfH3/69On+/n5bp+SY8vLyuLg45pJLfHz8xKZ5e3tbLzLLZDI8kqK3t7e5uTkuLo4Q4uLiEhUVNbFphJCioiJnZ+fVq1fbPCAX3bhxo6WlhWmaq6trRETExKYVFRVt2LChr6/v5MmTfX19bMTklvb29o6OjpiYGELIrFmz1Gr1xKZJJJLh4WHm2GAwuLu72zrllGJvX7E2EUVR3d3dCoWCealQKDQazbgajUYTEhLCHHt4eDg5OWk0Gp4/M1qj0dzetO7ubovFIhAIJlaOjo7m5eWtX7/etgE5R6vVOjs7W884CoWipaVlXM1ff/313nvv/fTTT1euXLF5QC7SarUikci6gqVQKK5fvz6u5tq1a8PDw4cPH/by8qqvry8uLg4PD7d5Ug7RarXu7u7WhdJJz2m5ublJSUmRkZEWi8XBweHQoUM2jzmV2P8nQoqiKIqyLiYLBAKz2TyuxmKx3L6vYdIavrm9JwKBgKZpi8UyaVlycrJSqczOzrZtQM65n3fRpk2b3n77bU9PT9tG4y7mNG19OWnTRkdHBQJBU1NTSUnJ+++///rrr9s2I+fcT9OOHz+u0WheeumlxMTE69evV1VV2TbjFGP/nwiFQqGHh4dOp3viiScIIT09PXK5fFyNTCazPqFXr9ePjIxMrOGb23vS09Pj7u4+cbctRVEpKSk3b9784YcfJv2wyCtSqXRkZOTWrVuurq6EkJ6eHplMdntBV1dXaWmpSCT6+eefu7u7h4eHMzMzc3JyxpXxilQq1ev1RqPR2dmZTNY0QohcLg8LC2NO/SqVatOmTWazWSi0/3PXnUil0v7+fusVmp6eHuvFG6t33nmnsLBw5cqVhJDZs2dv27YtKSmJhaxThP1/IiSELFu2rLq6mjmurq5WqVTMcV9fH/MpR6VSMfsaCCFVVVXz5s3D3+wqlcratKqqKmvTBgcHTSYTIYSm6ddee62zs/PIkSPMWYzn5HL5448/zjSNoqja2tply5YRQsxmM7OyJRaLv/76a7VaHR4eHhQU5OjoGB4ePnPmTJZzs0qpVPr6+tbU1BBCLBbLiRMnrE2zrtO/+OKL7e3tzHF7e7tUKuXzFCSEBAQEiMVi5ukCY2NjdXV1TNPGxsYGBgaYGoFAwPyeEkKMRiP+Tr0Hljfr2MT58+dFIlFOTs6bb74pkUi6urpommbeJefOnaNp2mg0zps3b+3atfn5+Z6enocOHWI7Mvu6urokEsmWLVtycnJEIhHTKJqmvby8vvvuO5qmP/roIwcHh6SkpIyMjIyMjKysLFbzcsJXX30llUrz8/MTEhKefPJJk8lE03RTUxMhZNytFA0NDdg1yigsLFQoFLt3746NjQ0MDDSbzTRN19fXOzo6MgVarVYul2/duvXDDz/08vL67LPPWM3LCQUFBb6+vnv27ImOjn722WeZd1dFRYVIJGIKdu3apVAo8vPz8/LypFLp7t27Wc3LdXx5+kRra+vhw4ednJySkpJ8fX0JITRNf/HFF7GxsczuhsHBwQMHDvT19S1fvjwsLIztvJzQ1dV18OBBk8mUkJAwf/585odFRUXPPPOMn5/fuXPnfvnlF2uxUChMTU1lKSmH1NXV1dbWPvbYY+vWrWO+mUGn0x07diw9Pf32su7u7urq6uTkZJZickttbW1dXZ2np2dKSgqzcaa7u/v48ePWd5RWq/3mm2+MRmN4ePiSJUtYDcsVlZWVDQ0NCoVi3bp1zHWFP//888SJE6+++ipTUFdX19jY6ODgoFKpnn/+eVbDch1fBiEAAMCkeLFGCAAAcCcYhAAAwGsYhAAAwGsYhAAAwGsYhAAAwGsYhAAAwGsYhAAAwGsYhAAAwGsYhAAAwGsYhAAAwGsYhAD2IzMz08fHp7Ozk3lpMpnCwsJCQkJGRkZYzQXAafiuUQD7MTQ0FBwc7Obm1tjY6OTktGXLln379p05c2bx4sVsRwPgLgxCALty4cKFpUuXZmdnv/DCC6tWrSooKNi8eTPboQA4DYMQwN58/PHHb7zxhqura3h4eElJCfNsdwC4EwxCAHtz8+ZNHx8fvV7/22+/LViwgO04AFyHzTIA9mbDhg3Tpk3z9vbOysoym81sxwHgOgxCALtSWFhYXFy8d+/e77//vqmp6d1332U7EQDX4dIogP1obW0NCQlJSUn55JNPCCEffPDBjh07Kisr1Wo129EAuAuDEMBODA8Ph4SECASC5uZmFxcXQghN06tXrz579uyvv/4qk8nYDgjAURiEAADAa1gjBAAAXsMgBAAAXsMgBAAAXsMgBAAAXsMgBAAAXsMgBAAAXsMgBAAAXvsvBYzI5L2dhdQAAAAASUVORK5CYII=", "text/html": [ "\n", "\n", "\n", " \n", " \n", " \n", "\n", "\n", "\n", " \n", " \n", " \n", "\n", "\n", "\n", " \n", " \n", " \n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n" ], "image/svg+xml": [ "\n", "\n", "\n", " \n", " \n", " \n", "\n", "\n", "\n", " \n", " \n", " \n", "\n", "\n", "\n", " \n", " \n", " \n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n" ] }, "metadata": {}, "execution_count": 6 } ], "cell_type": "code", "source": [ "rvecs = collect(r_vectors(basis))[:, 1, 1] # slice along the x axis\n", "x = [r[1] for r in rvecs] # only keep the x coordinate\n", "plot(x, scfres.ρ[1, :, 1, 1], label=\"\", xlabel=\"x\", ylabel=\"ρ\", marker=2)" ], "metadata": {}, "execution_count": 6 }, { "cell_type": "markdown", "source": [ "We can also perform various postprocessing steps:\n", "for instance compute a band structure" ], "metadata": {} }, { "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Computing bands along kpath:\n", " Γ -> X -> U and K -> Γ -> L -> W -> X\n", "\rDiagonalising Hamiltonian kblocks: 12%|██ | ETA: 0:00:01\u001b[K\rDiagonalising Hamiltonian kblocks: 91%|██████████████▌ | ETA: 0:00:00\u001b[K\rDiagonalising Hamiltonian kblocks: 100%|████████████████| Time: 0:00:01\u001b[K\n" ] }, { "output_type": "execute_result", "data": { "text/plain": "Plot{Plots.GRBackend() n=126}", "image/png": "iVBORw0KGgoAAAANSUhEUgAAAlgAAAGQCAIAAAD9V4nPAAAABmJLR0QA/wD/AP+gvaeTAAAgAElEQVR4nOydZ1xUVxPGZ5cFlCJtCx0EKxGlCKJgQbGDvUXFaEzQWFBjDOY1CilGjImisRuN2MWoEaPGYEERsWCv2LuiAkovu/u8HxZpLmzhbiHu/+cHuXvOnAEu97mnzAwLAOnQoUOHDh0fKmxNO6BDhw4dOnRoEp0Q6tChQ4eODxqdEOrQoUOHjg8anRDq0KFDh44PGp0Q6tChQ4eODxqdEOrQoUOHjg8anRDq0KFDh44PGp0Q6tChQ4eODxqdEOrQoUOHjg8anRDq0KFDh44PGo6mHagVmZmZb968cXFx0bQj0klLo8aNia3Kl42rV6lFCyKi7GzKySE7OyaNi8ViFovFYrHe/yg3l27cIB8fEotp3z4KCWFy3LLR2Sr92amA06fJx0e1v3HVAQCAEj/zJ0+eXbhQbG3tkJ5+qVu3FgYGBqpwT1Pk5dH589S+ffmV5ctpwgSid7foxYvUpAkZGdHff1PbtmRlpVp/kpMfHTzI7dvXqLAw39LyVfPmTqodrwKHDh1avNg4K8vv8WNWYWGCt3dXyfW3b0ksltH3xYvTenptSkooL+/FhAlvv/++qZyDFhcX//TTbnv74AcPLvftW8/Hx1NRtwFIfYhVpG7+yb5j//79s2fP1rQX0klJocBAysxU4RDff0+hoZSfT0Q0ejQtX86w/eLiYpFIJPWj69cpJIRu36asLPrmG/r6a2I2Zy2AfMk3Vqf48kvMmPHvtm3bcnNzNe2LwohEouLiYiU6AqIFC74cNKjtiRMH/mMqWFRE/frR1q3lV0aNohkzqLi4/BZduJC++46I6NQp8vSkI0dU5Ux+Pv3wA3Xvnvjzz5+1bn2id++w48cTVTWYNCIjI//+e2dyMuvRo+KXLydnZZHkn7ExWVjI+Pfw4Zi7d98+ekQZGcuvX18r/6AGBgYnTvw5YYLb8uWficUK359isVgoFMpuh7rMxo0bhw8frmkvpJCXhyZNsGOHCodYuhSNGuH5cwCIjYWbGwoKGB6ioKCgpKSkuk/Xr0fTpnjzBpmZCAhAaCiKixkbWiwW5+TkMGZOXQQFfc5mjzc0nOvi4pubm6tpdxSjpKSkgPF7qC4jFGLQIAwYAKGw9EpiIths7N4NVLhFX72CjQ3OnweAw4fh4IDwcBQWMumJSIS4ODRsiMGDcfcutm7dERo6ZetWVT5fpPH8+XMiAVEPoiaOjl4K9Z05M4rNdmSz25qYuKrzT1vyeiezmU4IVcKECfjkExXa37gRDg64fx8AnjwBn4/UVOZHqVkIAUyciO7dIRQiLw+9eyM4GHl5zAxdF4VQKBTy+b5EIIKZ2cyDBw9q2iPF0AlhRUQiDB+Obt3KJa2kBGZmGDSo9MuKt+i6dfDxKdXLrCwMG4YWLXD5MjOeHD4MT0+0aYOkJGYMKs3Nm+DxSqZPn33gwAElur99+/bChQuMe1UzOiHUGAkJsLdHZqaq7O/ZA2trXL8OAGIxevTA3LkqGUimEJaUIDAQ334LAEIhxo5FmzZ49YqBoeuiEAKwtfUgyiAS6ev3PX78oqbdUQydEJYhFmP8eAQGVlplCQqCpSVEorI25beoWIwuXbBkSXnj2FhwuYiJqZUbN25g8GA4OSE2FmJxrUwxQr9+WLAAubm5Ym3wRj50QqgZ3ryBoyOUemGSiyNHKs3/li6FtzeTa5IVkSmEAF6/hqsrtm0DALEYX38NNzc8elTboeuoEC5dusrAwElPz9HV9ctGjXDtmqYdUgSdEJbx9dfw9UV2dvmVLVvAZiMlpfxKlVv01i3weHj8uLzBzZvw8kK/fsq8Gr56hfBw8HiIjmZ4lVVpTp6EszMKCnRCqH1ooRCOHInJk1Vl/MwZ8PlITCz98u5d8Pm4eVNVw8kjhAAuXgSXW67NMTFwdsaiRTs//njSypXrRGWv0IpQR4XQ2dmX6DZRLpcbuHDhcz4ff/2laZ/kRieEEr77Du7uyMgov/L2LerVw6RJlZq9f4t+9x16967UpqQEkZGwscH+/fKOnp+P6GhwuQgLw8uXyvivItq1w6ZNAHRCqH1omxD+9RdcXcH4A3zFivU2Nh48nqeZ2ca//y69KBKhfftKqzGMI6cQAti5E87O5X+3kybtYrEGEJ00MZnw/fe/KjF0XRRCoVBobt5GskfYoMHs/fv3p6TA1haRkVqxriUTnRAC+O03NGqEZ88qXfT0hK1t1Zbv36JFRXBzw65dVVseOSLXCRqxGHFxcHZGcDDu3FHSfxURF4dWrUqXhXVCqHVolRC+fAlbW5w6xbDZzMxMHs+XqIiosEEDn7dv30quz5uHzp1V+4SVXwgBzJyJ9u1RVAQAw4ZNJDpFBKKXXl49lRi6zgnhjRvo1Qv16/fT1/+ZaKu5uWdmZiaAp0/h64thwxg7RqQ6dEK4fj0cHUvPoJWxdCn09KQok9RbNCUFdnbIyqra+M0bfPxxTSdoDh+Glxd8fXH8uNLuq4riYjRujISE0i//k0JYt+MItYrx4+nTT6lNG4bNZmZmCoVORAZEhvXqOWRkZBDR9eu0aBGtW0ey4kTVx9y5ZGZGM2YQEbVv72lisoUoi8WK9fT00rRrqiUjg6ZMoY4dqVs3Sk/f8ssvpmPH3qtff4+BgQUR2drSsWNkYEABAfTokaZ91VE9u3bRN9/QwYPk7Fx+8fFjmjqVIiPJ1VUuI35+1Ls3vR/bbGZGW7bQjBnUpQstXkwlJcKkpKQrV64QUVoaDRlCY8bQlCl06lSlyH0tYcUKatyYgoI07YdKUYMmqw7tmRGuXQsPj9L5ELNs3izmcLrWq/c/U9OZvr49xGJxSQlat8YffzA/VkVsbFoS2RLZBAT0kLNLdjbc3LBmDUQi0axZ89zdg1q3nunnl6/Ebn+dmBEWFyMmBlwuwsOrTgJGjsScOZWuxMTA3p75BQMG+ZBnhAcPgs8vjQWsiLMzWrWS3qW6W/TNG9jbIzlZeq979+DnV2RqGmhuPtHSsq+X1/8EAkRHMx8EzBTZ2bC1rTSRrUMzwpUrf9fXd1yy5B+ZLXVCyACPH0MgwKVLDJsVixEZCVdXXLpUvHv37j179kgWKmfNQkgIw2NVYeXKlUStiIqIsojs5NekmzfB45Uv74jFGDIEo0cr7ID2C2F8PFxdERQk/WjokyfgcvHgQaWL+/dDIFD5G4zSfIBCWFxcXFJScuIEeDycOFH104gIGBhUe+azhlt02za4u1d7lvvgwcP16k2T7CUbGbVOT5d390EjfPMNxo6tdKUOCSGH40D0ZONGocyWOiGsLSIRAgPx888Mm83ORt++aN++6smx1FTY2JQmlFEdo0ePJhol+Vsl8ktKOi1/33/+gZ1d+TnynBy0bImlSxVzQJuF8Nw5dOgANzcZQTKRkRgxourFtDQ0a4awMMi996o+PjQhnDNnAZfrYW7eyth48dGjVT89fx5sNtatq7Z7zbdoSAiio6V/dPr0aSurEUQgKhAIWmmzqDx5AiurqtFQdUgI2WwXIgweLPvguk4Ia8uiRWjXrjwJEyPcvYuPPkJYWNWXyoICfPQR4uKYHEsq8fFPiQREPxF9SdTZykqxkN7IyFxT00ECgbe3d49nz549eABbWxw+rIAD2imET58iLAwCAWJiZP/G8/Ph6ChlnpGRga5d0b27lCMVmuWDEsL09HQuN4BITCRq0MAvq/IvQyQCj4fAwJos1HyLPnwIHq/aw58jR07mcn14vJaxsduV8V5djB5dmi6jInVFCB88AIs1imjczz/LPoOrE8JaceMGuFzcusWkzaQk2Nhg8WIpH4WHS5lkME5KCvT10b379aZNm/v4+Jw48cbODgIBWrfGsWNyWfj223ls9jIisFgHBg0KA5CUBGtr3L4trw/aJoR5eYiOhpUVwsPx7tyubGJj4ecn5QVCKEREBBo3Lk0PpCV8UEJ47949S8v+kjUPPr/344qR8MDQoTAyknHQV+Yt+ssvCAys9vWxoKBAyOzrM9NcugRrayl3e50Qwjt3YGwMd3fx7NlzLl++KrO9TgiVp6QEvr5YuZJJm6tXQyDAoUOVLl65ciU4eHS7dp8KBNcqxvmqglOnoK+PwYOBCuETGRno1g1ubnB0RFAQrsq6r0aOnEKURASiJ76+pfuZy5aheXN5VUR7hLBidNe9ewr39fHB1q3SP920CXw+4uNr7yMzfFBCeP682MBgoInJZ1ZWn3bv/nHFj44cAZstOzmUzFtUKIS3d2kQel2ke3fpOxraL4SXLqFePfj5Abo4QjXw3Xfo2pWxSL6yWcKNG5WuFxQU2Np6Ep0mSrGy8ixSxcnUd5w5A3199OpVPnRZHKFQiMhIODjgyy9L016kp1drJyXlFJfro6f3G5vdeeHCP8uujxuHvn0hT6oZDQphSUnJ48ePJd/46dPw90fr1srnOz55Eg4O1c4tTp6EnV21m0lq5sMRwvPnYWODnTvFycnJKSkpFR/rBQUwMcGwYbKNyHOLXrwIa2vtShAjJ0eOwMVF+jF4LRdCyYJW2bK2TghVy/nz4PEYSKopITMTQUHS941u3LjB50u21iEQfHyL2XXYCly5AkNDdO9efuX9gPo9e8DnY+lSRETAygqRkdUe+759+/Yff/zx9deXfHzKdzqLi9GpE2bPlu2MpoTwzp07jo7eAkFfOzvv/v3v2ttj1Sq5lLsGBg+uKSv6kyfw8cHHHyM/v1aj1J4PRAhTU2FtXW3qu06dYGXF5Lva1KkYM0ZBFzWNSARvb+zcKf1TbRbCY8egr49+/cqv6IRQVdy+fTs19bK7O7ZsYcbgrVs1nSQsLCzkcr2IjhEdtbf3kueXqgTXrsHQEF27VrooNbPMjRul3l6/jsGD4eBQk1SIxejdG1FR5VckSbqrWy2s0FEzQjhkyBdER4lAdNTdfQIj4nTvHrjcqlm7KlJQgNBQeHnh9OnHv/yyOC5uh3IJWmvJhyCEZ8/C2hp79kj/dPNmsNk4c0YuU3Leonl5aNiw6maHlhMbizZtql3r0loh3LcPenpV69/phFAljB07nccLNjIaxucPZeRuOHgQAgHWrq22QVYWbGxuduw4bsiQL1Q0HZSoYFBQ1evVpVh7+xZ9+6JDB7x4gZSU0sXDslTgVUhPh50dTp4svyJJ0n32bE0uaUoI+/X7jCiFCEQp/fp9xpTZiAh8+mlNDcRizJr1XE/Pg8Vaa2IybdiwL5gaWn7+80J44gS4XOzdK/3TrCwYGmLKFHmtyX+L7t+Pxo21N2S+CgUFcHKqaS9AO4VQUh4kLKzqdZ0QMk92djaP106ySmllNfzKlSu1NLhqFayt8X4MU0VGjkR4eC3HqYkbN1CvHjp3lvJRDblGxWJER8PREadPlx4ncXFBUBCk/kj274era6VjMrt2wdm5pl1GTQnh2bOX9fQ8zczGWVt71P73W4YkPUfNxZO3bNmir/+L5O6ytm7N1NDy898WwqQkCATlCTPfp2VLNGyogEGFbtFBg+TaEdAGoqMxcGBNDbRQCDdsAJuNGTOkfKQTQubJy8uzsGhNJCICl9snLS1NaVMlJZg4Ee7uVTP8VmHPHri4MF/Ooox792BkhI4dpX8qM+n23r3l09miIsTEgMdDWBhevMDp06ejouYferckFBZWNb/M7Nnw9682H7+mhPCHH9C/f/bZs2cZH331anTqVFODkydPWloOIiohSrO3D2B2dHn4DwvhsWNSDmNXZP586OkpdipYoVv0+XPw+XWgPmVmJvh8GVE92iaES5aAza52G14nhMzz6hWsrH41NfXl8QLGjPlSaTuvXyMwEL16yYglePUKtrYqzEb/4AGMjODjU20DeapPpKWhefPy2P+MDEREoEGDw/XrdyTabm7ef/nydQDy8tCsWaWtwZqzr2lECF+8AJerQLCjQohE8PLC7t01tZkz52d7e29X1y6WlpdrXidQBf9VIUxMhECAI0eqbfDgATgc/PKLYmYVvUWXL4efX20PXqmaqVOrll18H60Swh9+AJuNBQuqbaATQobJz4efH77/HtnZ2ZIKOwpRUFAwZMh4BwfvwMAxTk454eGy/ySGDMFXXynprUwePYKxMby9a2ojZxmm7Gz074+AgPLEbyEhk4lOEoHoRevWpbVKz50Dn4+HD8s75uTA3V16rJJGhHDUKPzvfyq0LzmSLk8K8qQk8PnKx2wox39SCP/5B1yujN0HR0d4eSlsWdFbVCSCvz9WrVJ4ILUhOdX14oWMZtojhHPmgM3G6tU1tdEJIZMIhejXr2ryWYWYNWueoeEvRGCxVnft+o3M9rt3o1kzVR2pl6igp6eMZvLXI5RsGTo4lFZXmDcvpl69uURgsbYMGFC+w/nTT+jQoVJysvv3pWdfU78QnjsHGxsFssYoR0iIvDMPySmqmrcVmeW/J4SSLOcpKdI/nTNnrr6+E5vtxuHsVSJPhRK36NWr4HLx9KnCY6mHoUPx44+ym2mJEE6eDDYb69fLaKYTQiaZNAk9etQqUXL37mOJLhKB6H5AQI2b0e9q/J5WINO1Ajx5AhMTeHjIbqlQYV4Af/8NgQBr1qCoqGjYsC/s7b2bNx/u4pJZ9pcvNUG5JPtalayMahZCsRjt26ujLsSdO3K9dEvYvRs2NrLz+DDFf0wI9+2DQFBt3atnz56x2Q2JMolesViuGYoroXK36KxZGDpU0U7q4PRp2NkhN1d2S20QwvHjoadXbaRjRXRCyBjz5sHdHW/eKNldLEZsLMzN9xoadiPabW7e948/NtfcZeBAVa3RpafDzAzNmsm1V6GoEAK4das0XfjJk2fnzv3lyJEjixahRQu8fl3a4PFj8PlVQ7WWLq2afU3NQrhlCzw91bR/M3UqJkyQt/GOHbC1rZpsSEX8l4RQcoyrhlfJAwcOsNkhkjO6bHa3Y3Jm0a2AcrdoYSGaNdOivHpldO5cUxBXRTQuhMOGgcPBv//K1VgnhMywbRscHFA5Ja8CXLiAtm0REIDLl5GUlDRz5g+HDsmowrBxI9zcVBJ19OoVzM3RtKm8T3wlhBBAdjbatz/C4XQg2mJu3nflyvUzZ8LDozxpzo4daN68ataxKtnX1CmE+flwdpY3n3jtycyEQFCp0mnNrF8PR0eF05wqwX9GCHfsAJ8vI1C1oKCIxWpJtJJomYGBvRL3udK3aGIinJxEiYmnVZclSlH27EHz5vKueGlWCHv1AodTbdTy+2i7EGZmZtZw5EQsFr98+VLmeoWqhfDYMfD5SlbcffMG4eHg8RATo8BU4+lT2X/DyvHqFSws0KSJAgu8ygkhgOHDqx6WmTYN7dqVL7yMHFl1ViTJvlZW1V2dQvjdd+perVqyBF26KND+99/h7Fy1zC/jaEQIhULhrFnzvL17ffllJCN5dLdvh62t7L/ZDh1Qv/6roKAhISFDHimVKVHpW7S4uNjKqmv9+p9xuT2mT/9OCQvMIhTio4+wf7+87dUvhC9evPjpp58OHDjQsycMDBTbONdeISwqKhoyZIiVlRWXyx08ePD7d39ycjKXy7WxsbG0tGzevHlq9d+3SoXw2jUZsUc1EB8PR0eEhlZb3ro6evWqlJCs9rx69So+Pv7hw7cWFmjUSLFtTqWFcO7cRfXq/UQEFmvr8OGTAYjFGDsWXbuWnpnMyUHjxlUXiCTZ17ZtA9QohJJS8jVHczJOSYlijx4AixahSRPVFmTWiBAuWLDU2Hg60Yv69aNmzqw+Jat8bNwIW1vZu6pTp4LDkZ78QX6UvkWPHj1qbh5OBCIxn++t8WJMq1ZJz6dRHWoWwrS0NA7HlsWaStRRT+9XRbcJtFcIV6xY4eXlVVBQUFhY6OPjs2zZsioNMjIyXr58CUAkEkVERLRuXW2WDdUJ4bNncHbGhg0Kd7x9G927o1WrSknF5GTtWnh4VC3GWxu2bdvGZjuw2f2I3OzsrioqakoLYVFR0dCh4+3tvQWCjz/5pHTeLxRiyBD0718qxidOwNq66pNdkn0tNVV9QjhypJTSo2pg3z40a6bY7zoqCu7u5butjKMRIezV61OiS0QgemRh0f/rr7FpEy5elF73oGbWr4ednewSjytXgs3Gjh3K+VuO0rdoSkqKldVoIhAV8XgarlCfmys759F7XdQqhMOGjSJaKvlx6ek5K9pde4XQ399/5bsifr///nvbtm1raBwfH+/q6lrdpyoSwuxseHoqXBwnPx+RkbC0RHS0MgXrnzxRbOtIHqysWhJdIALRwUaN2ivaXWkhLCMnB25u5ZvwRUXo1QujRpWuFc+eje7dq+b23bkTzs6YP3/55MmTnzx5UpvRZZKaqo6Qieqort5bDXzzDTw8oHgUq1yoXwj37gWXu5HDGUl0zth4/Oef/xYTg9BQeHvDxAQuLggORmQk4uJw9ar0HNBisfjhw4fZ2dlr18LWVnbqlqQksNnlK/C1QWkhFIvFAwd+xuN1qF/fo2XLDZo9gBkZiZEjFeuiZiGcNCmcaCYRiJ4bGDgp2l17hdDW1vbouwDX48ePW1tbv9+mqKho1apVkZGRnp6eu3btqs7Uxo0b+/Xrd/cdyq31V6G4GN274/PPFeu1dy8aNkRwsJKFmcRi9OiBn35Spm8NcLmtiM4TgehA48YKZ+2qvRACSEsDn1/+ypmfjw4dMHEiAJSUwM9PihgIBF2IQlisORyO3QOV7YyJxQgIQGysiszL5vp18PkKz/C++gp+fsjOZt4fdQrh48cIDUXjxvjnH/zxx6b+/cNWrlxXseZGcTGuXkVcHCIjERwMFxeYmcHbG6GhiIlBQgJevkRhYaGPT3c+v4+pqTeX+6fMlEBPnqBevdKi07WnlosWWVlZublFbdti0SJm/FGC9HRwubh7V7FeahbCgoICNrsni9VKT89u6dIVinaXUwg5pHays7Pr168v+b+xsfHbt2/fbwPg3r17mZmZeXl5eXl5NVg7efJkly5dJP93cHDYv39/Ld2bNs2QiPXzz4W5uXK1f/aMFRVlePo0e+HCoqAgERHJ2bEia9fqv3ihP358vhJ9a6Bbt++3bBnDYjmx2anLl8fmKmi9sLCQw+FwOLW6SWxtackSzoABBsePF1hZgYi2bWMFB9ePiBDOnl28Zg07MLB+69YFH30kLuvy6tUdovsASyg0WLBgQXR0dG0cqI64OE5+vkG/fgz/zOXHwYFCQgyjomjevCL5e82ZQ9OmGXbvzv7rr0IjIzDoj/AdDNp8n6IiWrLEYNky/fDwkpiYYgMDIuo7aFBfIsrPz6/Y0smJnJyoZ8/SL7OyWFeusK9fZ1+8yN60iX3jBrt+/fisrDYi0XdEeSYmnaytu9fwqywuJk9P40aNsG4dM79xAFUcVggOhwMUr11bEhho1KJFoZ+fiAGfFOSbbwxHjACfX6zQDyQ/P18sFrNYLJX5VYnVq/XZ7P3Xrz+1sTEjIkUfYmKx2NDQUHY7xRW6tri4uPzzzz+S/x86dMjZuaZl35MnT9avXz+vmgrfjC+NfvcdWreWK6oUQHFxaZrpyEi58mZVx/374HKZj5u+dg16evj++7f79u2r7gdYM4zMCCVMn46uXctXjF+9wkcflS4+r12LFi0qhYtwOI5EL4hAFLpw4UJGHKhCXh4cHXHihCpsK8DLl+DxoOgp+ionj5hCDTPCf/9Fkybo379Spj3lEIkQE7PNwOAHIhAV2NvLqNfh7g4ul8moJKa2sf/+G46OCp+qqz03b4LHgxIpddQ8IzQ1lZ3+tAa0d2m0X79+P/zwg+T/8+bN69OnTw2Nnz59SkTVxVEwK4SbN8PFRUbWj61bd9nZeVlbe0ycuKJFC3Tpgps3azWoSIROnRRO+CuTkhJYWVVbVkJOGBTCkhJ07FjpQOzTp3B1xfLlADB0aKWsqosWLdXTc2CznevVGzd1KiPjVyUyEh9/rBLLirJgAfr2VbiXUIhhw9CvH5NHq1QqhE+fIjQUrq7Yt48xm/n5+e7unXi8ETye35o1G2to2b8/DA2V3LaoDgbPc82YgV691J2Pu08f/PqrMh3VKYQREahXr1YpvbRXCA8ePMjj8Y4ePZqYmCgQCA4cOCC5HhAQcOrUKQDbt2/fsWPH5cuXjxw50q1bt65V6qZXgEEhPHIENjYyVK24uJjHa0mUS1TMZndcs0bZMPsKLF6Mtm2VOVxTM127wsxMmaN3FWFQCAG8eAEHh0rPwYcP4eSEtWuRlQUnp0qRKpKnTFYWvL2ZOdpQEUnIhKrD8uSkqAiNG8ubJqMixcUICcGgQbV6TFRERUJYUoKYGFhZISKC+TQRIpHoypUrkkPm1SEpUFBD9QnlYFAIS0rQvj3mzWPEmAzy8vJ+/XXp8OHf29vfV25FQW1CWFQEQ8Nq6yvJifYKIYDY2NiAgICAgID1FXKmjhw58tKlSwASEhIGDhzo7e0dGBgYGRn5pvrkZkwJ4dWrsLaWkew/MxOLFmUYGAS9K8z7eQ0BjnJy9y74fIVXxmSybh3Y7GpzDcsPs0IIICUFAkGlzfm0NNjaIi4OiYlwcCg/OVL2lHn1Ci1ayJULWH6GD2deXGvDzp1o1UqZl6GiIvTsiU8+YWYyoQohPHYMkoUT9WSJe5/4eNkFCpSD2Qif589hZ6dk1LJCBAT0NTT8lSjOwsLrtVKxOGoTwk8+gZlZbY1otRAyBSNC+PQpnJ2rjSsqLER8PEJDYWGB4GC4ufU3Nf3KyCjqo486yPPzrQGRCO3bY8mS2tiQgqS42jey61vIhnEhBBATg1atKuVXu3gRAgH278f06ejfv/RixadMejqaN1c4mqU6UlJgby/vNrDaCArCmjXKdMzLQ8eOGDtWenSBQjArhM+fIzQU9vaaPJd76RI4HHzxhUqMMx7qevAg7O3lTciuHEKhkM9vI3mVNzWN/Pvvv5Uwoh4hzMoCh8PAG4xOCOXi7Vu0aoX3D2SIREhKQng4+Hz4+yMmpu+oeJEAACAASURBVHSyIhQK4+Pjd+zYUftHxoIFCAhgfmPA1laZ4mpSUYUQAhg1qmroUkoK+HwcPgwPD6xbB7z3lHnyBK6uCkfdvY9YDH9/bKxpO0kzXLigfETj27fw9cWUKbX1gSkhFIkQGwtra4SHQ+01JcvJyICJiWI5UxRCFTkf5sxBYCDzGyVliMWwsGhDdIMom8sNTEtLU8KIeoSwd2/Y2TFgRyeEMsjLyysuRteuGD++0vWrVxEZCRcXuLkhMlJVyY5v3ACXW7X2UO0ZPBhGRow9fVQkhLm5aNGi6rveoUPg8xEXBx4PaWlSnjKPHsHFBSsUjiOqxIYN8PbW0irhY8YoP49/8wbu7o94vE58vnfbtsHZSoUZKi2EYrE4Lu7PqKj5169fT02Fry86dVJf9SipiERwcoKDgwp/16oQQpEIXbsiMpJZq6WUlGDsWLRsecXTs1fjxv5//LFFOTtqEMIHD8BmQ6n5alV0Qlgtb968adWqM5/f3tjYOzDwuuT968kTxMTA3x/29ggPx7lzzHtbRkkJfH3xLrsOY+zeDTYb784eMYCKhBDArVtSijH99RdsbDBrFlq3RlGRlKfMnTuwt1dyCRFAfj4cHdVd9l1+XryApeXz775bGR8fr8SDpmfPT4kOE4HDWfO//ylzwEBpIZw8eZap6RdE2+rXby0Q3IiNZWCdtpb4+6NBA9XmDFJRFsD0dNjb4118GWPk5qJnT/Tty0CtbzUIoZ8fmjVjxpROCKtlxowfWKzVRCC66OUVEhuL4GBYWSE0FAkJ6pguzJ2LLl0Yflikp8PQsOrstpaoTggB7NkDJ6eq4VMbN8LBAR4eBw0MmnA4DgMGfFKl161bsLfHZhn1HKUzezZGjFDSWzWQnp5uaurBYq0wNZ04alS4ot3btOlD9IgIRInGxlOHDsXChUhOVuCgphJC+Po1/v0X5uY+RGIiEG2ZM4fpSCDFmTIFHI7sdGu1RHXpcBMTYWMDBtMLvngBb298+ikzB4xVLYQXLoDFqraisqLohLCcZ8/wzz+YPx8jRsDdHWz2N0R7JcnrWKxmAwZg1y6GY5Nr4No18PkMBBRXwdUVTZowbFOlQgjg66/RpUvVHZFly8BiORFdJCpgsTpsfG9D78oV2Nig+rx70nn4EFZWDEeSMcvWrVsNDX+RHGSwsPBWtBD0zp17rKw61Ku3gMfz+uefc3FxCA+Hvz+MjeHmhtBQrFqFq1dres+TRwifP8e+ffjhB/TvDycnNGiAjh1hY9OHKJlIZGr6+e7dfynmN9MsWwY2W67a5bVEpXnh586Fnx8zQaJ376JxY0REMGBKgqqFsHlz+PkxZu1DEUJ7e7eUyoECZVkK//c/dOwICwvUrw87O7i6wtERDRqAw5lE1JxoJpGvo6OMhBRMUVhY+MknUxwdW1tafrFiBcNnFqdMgYEBnj1j1qrKhVAoRLdumD276nUiR4keEP3cr1+/9zuePw+BQLEthI8/VtXWC1OcOHHCwmIokYjonpGRv6kp2rVDVBSSk+V9kb958+bWrVurZNwtLkZqKiTJrN3cYGYGf39ERCA+HhWj7y5cuMDlNmrQoOH27dsrdn/6FPHxlRJ++vsjPByxseWy+vDhw3bt+jg4eGu8ut6xY2Cz1fSLro0QvnnzJjY2ds+ePaJqXkzEYvTty4B6nT0LGxusWlVbOxVRqRAePAg2m8nDEx+KELJY29hsx99+ezRhAvz8YGMDDgf168PICBwOrK3h74+wMERHIzYWCQm4exdXrlzl8dwNDT+1sPDbsqXW5VjkY86c+fXq/UQkZrOXTp48i0HLhw6BzUZcHIMmS1G1EAJIT4eDQ9XpHZEdURTRn0SePj7zpXaUaKGcaUpOntTGkIn3iYj40dbWq3nzjufOnS8oQEICIiLg7Q1TUwQFITpasXI5UsnKQkJCqbBZWsLGBoMHIyYGRNZEvxFtJLJbufKRpIFAAAuLSsqn8c2/GpDk1FZbjWWlhTAnJ8fZubWBwQJT0ykhIZ9U1ywzEw0bYvdu5T1MSIC1NfbuVd6CVFQqhPb26NWLSYMfihASgWiWvv52Jyd06oRJk7BhA06elFHC9NmzZzt27LhZy/RoihASUlZ37UFAwECmzL59CyMjVWULU4MQAjh1CgJBpXdAY2MXonCi6UQhfH5eeLj0Bb0TJ8DnIzFRhn2xGL6+2LSJSZ/VzIsX2LQJo0bBxgYuLhg3Djt3IiurVjaFQty7h3XrMHo0GjdOJfJ6NwufZGu79bvvsHcvnj5l6BtQPXl5sLSEp6f6RlRaCA8cOGhsLKkrBD7fr4bH9OnTEAiUPLi+YQMEApVk01WdEG7YAD09htOufihCyGJlslgeRxhPoMQ0MTG7WaxeRHvNzQeuWaN4wd9qaNWKmWgbqahHCAEsXYqWLcuj7C9duuTk5Gti0iwkZAmPh1atMGiQ9NNuR4+Cz0dyck3GY2Ph56fVUxmFuHIFv/6KHj1gago/P8yejaQklJTg7du3Awd+7uzsExY2o8pv7fVrnD2L7dsRHY1x49C1Kxo1gqEh7O3RoQM++QSzZxcQ2RA9JMokampkdLFHD6xbp6rCh6rAzQ08Xm1zCiqEEkJ48SK++gp8/hUOJ5hISPTK0NCz5lR/CxfCx0fhEwwxMWjYsLZpkKtDdUJoaYlRoxi2+aEIob5+07CwaZp2RAYZGWjcGNOmHf7yy8h9+xiLb4iMBIfDfDBiGWoTQgCjR1c60ln2lDl0CDweOnRAmzbSk24cPCglEqOMnBzY2zN2Ak2rKCzEoUOIiICXF8zM4Oj4JZu9lkhkaPjtgAHLJGl6PDzQoAHMzeHlhYEDMWMGli/HgQNIS6v6bA0N/ZTIlsjW17dzQUFpNiXJjuCqVVD05I56EIlEjRq1Y7Od2Gw3ff0kNZ+Ekl8IJaFZXl5wcEBEBG7eRGTkLzY2Hs7ObaZNS+RyERNT7YuaWIyBAxEu9yFisRhffYUWLfCYgUTI0lGREP70EwwMVJKN9oMQwoEDGVtmVBHFxejcGTNnMmz2zBmw2cwHI1ZEnUKYnw9Pz/J4+YpPmcuX4eiI/v3h4IALF6T03bsXAoH0j779VuEC3HWRly/RsGGvd7WrUho1mvTzz9ixA+fOKTCru3PnzsWLFyteycnB5s3o2xdmZujXD1u2aDJTzPvMnTuXxQolAtEjDsdVzaPLFML8fMTFVQrNkioft2+jfXt07VrtSfLsbDRtKlfIUGEhhg5FYKBqX1xUIYQiEYyMMGMGs1Yllj8MIWzQoEvv3qHqrI+lKOPGoX9/hsMTi4pgYcHwrvL7qFMIAdy+DR6vdFejylPm3j00bYq+fcHjSd/837ULNja4cqXSxUePYGXFfKSKdvLrr8vNzEYQ7bOwCNq//6Ci3Tt1GshmN2WzPZycfN//9M0bbNiA3r1hZoaBAxEXB6UKXDJJQQFatPgf0SIiEIn19BwVtbB+/TZ7e28Hh9ZbtigTbFGdEAqFSEhAaCgsLREcjNhY2T8roRDR0eDxsGqVdLG8dAk8nozE5VlZ6NgRAwcyP6mqgiqEcPJkGBurJIb7QxFCI6MMS8tPT548qWlfpLNgAby8mD+vGBAALpex+jvVoWYhBLB3Lxwd8fKllKdMRgb8/dGzJ2xssGyZlL47dsDWtlIY9ZAh+P57FXusTezc+deUKbOPHJF1fOg9srKy2GwXIhERWKxekZGHTp3CrVtSqrZmZUGSgKJBg9KnfG4unj171qpVFwuLllOmqOCV/j0ePUKfPuBwYGz8jMVyJYpmsUK8vLopZCQrK4vLbU2UR5RrZeWVna3wn+j7t+jVq4iIgEAAb2/ExKDG2lBSuHYNvr7o1k36q9vy5XB3rzYvzLNn8PBAdcfKmIVxIczLg74+8zVZJXwoQsjhFLLZE6ZOPabxV9T32bMH9vbML9YvXgw9PVy6xLDZ91G/EAL45ht07oySEimv23l5CA5GQACaNkV4uJTcxLGxcHDAnTvIy8s7dqzYwUHzExdtpqQESUmYMwc+Pq+ImkrOMbJYQzw8kn190agRzM3BYsHKCo0bw88PvXph5EiEh+O77zB/PsaPR5s2MDaGgcFAopVEN1ms9itqmQ22Rg4dgpcXWCw4OUFSwO3OnTsTJ05cJvXNqDKPHuHQISxbhvBwdO8Oe/s0FmuY5FtmswdwOPdtbODlhd69MXYsZs/G0qXYvRsnT+LhQylnVVavXquv78hmO7ZoEfjoEaKj0bQpnJwQEYHbt5X/BktKapoahoYiLExKr+vX4eSkvjBZxoVw6FBwuQzaq4ScQsgCQHWWTZs2TZy4p2HDInf33QkJeqNH05QpZGOjabeIiOjaNerSheLjydeXSbM3b1KLFjR3LkVEMGlWKoWFhRwOh8PhqHykCojF1KLFirS0eUSsRo0c09KSKn4qEtHEiXT2LFlakp4ebd9OZmaVuq9bR+HhXxgaXszOzv/ss29WrBimTufrBPfu0aFDdOgQJSSQpSUFB1NICE2bFnT9ei6RmYXFk5cvr7DZbEljsZgyMykjgzIzy/9V/DI9nS5ebEZ0k4iItgsER7dvX9mxI8M+L1pE8+fTq1fUti0tXUoeHqXX586NWbNmi0DA37w5plGjRpKLOTl06xalpZX+u3WLbt2iBg2oaVNq0oSaNKFmzahxY/GAAYF37/ZgsURNmx49e/bQy5esFy/o2TNKT6enT+nlS3ryhF6+pKdPKT2dTEzI2ppsbMjGhqytaeFCJ5HoJJEt0ShT00/GjAkaMYKxv/SrV2nMGOJyac0asrcvv56bS76+NHMmjRpVfvH0aRowgBYsoOHDmRldJnl5eUZGRiwWixFrL1+SrS1t2KAq/8VisUgk0tfXl9FOVUKsFjZu3NinTx/J68nduwgPh4UFQkNVdW5Yfp4/h5MTtm5l2KxIBIEA/v4Mm60OjcwIAXA4jpJzHyzW0IXvl8gCoqPh4lKaMK/K6fOUlBRj41FEICoQCFpp8+axOnn9GnFxCAuDszMEAgwejFWrqoYJJiYm7lI0cx0AwNbWk+g3omssVicbm716etDXR/PmmDattjntcnIwbRpMTKCvj/79q0aYJSUlWVgMIiomumBrGzRuHAIDYWsLIyN4eGDIEHz7LTZvxtmz0rNv5+fnb9iwcdOmzfJkWH39GlevIiEBGzZg/nxJFsDS5Efjx0+s1TcpjZIS/Pgj+PzSqmRlXLkCPh/XrkEoFALYswc8HpN59uWB2Rlh585wdmbKmBQ+lKXRKrlGX75EZCR4PAQHM1CiXTkKCuDnh7nKFACQQe/eMDFR33KfpoSQzXZ6l8R57pRq6uytWwc7O0yfDju7SrEThw8fNjefJjlAwed7CVVX203LEIvF27f/OXnyrISE0jLnJSVITUV0NIKCYGZWnpumuoeY0tUnXrx44e3dzcqq1YwZpSmTEhMxciQEAhDByAgdO2LjRsW2r27cQM+e0NNDgwb49lvpfdeuXaent1QiSMbG3suWISEBDx+qI2a0efMOLNYQovl6enbKVfWTh6tX0bo1unev9D6xfPkLQ8NOXK6vnV1HG5vnKq2TIxUGhTAtDWw2VBoE/oEKoYTcXKxahcaN4e+P+Hi1BlOLxRg+HEOHMjnohg0bzc0/MjJqzWanqrOKkKaE0N29M4vVl+g7Is+tW6tNwh8fDx4P338PHg9b3tVWKyoqatWqs5XVOC43eNq0SPU4rA0sWLCsQYNRRAfMzLpPmnRg8GBYWMDbGxERSEiQKyKb2Qr1EjIyMH8+/PxgZAQ2G/b2GDkS58+/36z8WM7Bg2jVCiwWXFxkJBi7e/ehvr4nUZyx8dcffzyBWc9l8vPPP48dO/aO6sJ4AbzbNeTzy/OFfvrplyzWX0Qg2jtgQK3LMSsOg0Lo6YmWLRmxVC0ftBBKEIkQHw9fX7i7IzaWmVTuMpk9G+3aMVnLIiMjg812ILpCdInFcn2r0hprldGUEAJYvnz5pEmTdu58XnMStdOnYWODb7+FgwPmzSu9KBQKk5OTr1+/rh5XtQRv715E6WVxhDt2SDnzWTOqEMKKXLqEsDDY24PFgqEh/Pwwfz7OnLlqYODAZrfQ17f/4otUgQBsNjp2hDyzrOnT0anTrVmzfoqN3aT+G1Wl1SeqcOUKWrdGjx54/BjBwWOILhCB6HKPHkwnYpEDpoQwORkslsoP/emEsJykJAQHw8YGkZGlGRqLioquXr2qXCHvGoiLQ8OGSE9n0uaBAwfY7JB3x9t6/vvvv0xarxENCmHZUyYhAVwujh2rtuWdO2jcGF9+CS8vjBmj1iRb2sPr13Bxmc5irSIqMTL6ZuHC5UoYUbUQllFQgKVLS0tEEX1JtIsIREdZrJHjx8sbsy+pZ/n6tYp9rR51CiEqTA2/+irJ1NSdw/mfqWnLY8eOq82BMpgSQldXdOhQezMykFMI2So5qaNlBATQ3r104ADdu0eurvTZZ89cXNp17vxz48adTpw4ydQoqak0eTLt2UN8PlMmiYgCAgKANKIDRPvY7Cvt27dn0rrWExREW7bQkCF0+rT0Bq6udPw4JSaSuztlZlLnzvTsWUliYuLly5fV66nG2L+fPDyoY8fIAQMuuLgEjBlD4eFhmnaqJurVo4kT6cQJys0la+vTRAIiIhKYm19asYJMTGRbuHOHPv+c/vyTrKxU7KvWwOFQRAQdPkxxcbcLCpoKheZA0ytXbmnaLyX56y+6f5+2bNG0H2WoXJFViRIV6u/dg6dnFFEcEYju+foGM+LJgwewtVXJ8a3168Fm33Zw6NKsWfvkmtNLM402zAglHDwIgQCnT1fbPicHPXqgXz/MnFloYNDRzCzcyqrfxIn/U4evmuPNG4SFwcUFR48yYE1tM8KKbNmyhc12YrGms9mNFi9eKk8XSTa+1atV7ZoM1DwjLKNLlxFEN4lAdLtz52Hqd4CRGaFAgP79GXFHBroZoXQaNqSgIBaLJSYiIlFaGmv/fhKLa2UzJ4dCQmjmTOrRgwkXK3D/Pn3+Oc2Y0ejRo0M3bhxv164dwwPUEbp1o/XrqU8fOntWegMTE9qzh4yNac+eE3p6rd++XZyRsSsu7qBIJFKvp+rjwAFq0YIMDenKFerUSdPeKMvHH3987do/339vfuHCzvDwifJ0mTCB3Nzo889V7ZqW0rGjZ/36G4kyWawNlpZemnZHGZYvp4wMWr9e035URA2arDqUmBECePHiRcOGPnz+MGtrz7lzT7Vti0aNEB2tZIE3kQghIRg7Vpm+Mi1bW8NXSupHNaE9M0IJ+/dDIKipOK1YjLFjU1isPkQgymexWn7+ufjsWdW6qn7evsXYsWjYkOFz5xqZESrKihVwd9eKhEGamhEWFxdPnx7VsmXXTz6JcnQsXirXLJpJaj8jbNAA48cz5Y4MdIdlaqKkpOTWrVv57zL3paaWJskND8f9+4qZmjoVXbuqJO1n9+4wMdFkvn9tE0IAu3aBz0cNsVOJiYl6en5EPkQtTU29580Tu7rCzQ3R0QqfotROEhLg5ISwMDB90qsOCOGFC+DzZeSeVhuaEsKKPHiAZs0QHq7WCLFaCuGsWahXT+V5ksvQLY3WBIfDady4cf369SVfenvThg109SpZWJCvL4WE0KFDctlZt47++Yfi4ojxNGRr1lBCAh08KNfZgQ+H/v1p+XIKDqarV6U30NPTMzNrSXSc6DSRcMYM3L5Nq1bRvXvUqBENGSLvb1YLyc+nmTNp9GhaupRWrSJTU007pF6ysmjgQFq+nJo107QrWoOTE508SefP0yefUEmJpr2RA6GQFi6kmTOZf2DWkg9UCKViY0NRUfToEQ0eTFOnkqcnrV5NBQXVtj9+nL79luLjydycYU+uX6cvvqAffqAPdU+wJgYOpCVLqHt3unZNyqft2rXz8cnh8/tbWnbgcicOGMDOzKSAAFq1iu7epaAgmjaNmjen+fPp9Wu1u14LTp4kDw+6d48uX6bgYE17o3bEYhoxggYPpoEDNe2KlmFhQYcOUUEB9epF2dma9qZ6duz409DQ2cCgSVHR+Nmza3coQxWoYXKqOpReGpUHSfShQICICDx5L73JzZvg81WSHKikBFZW6oiwkYkWLo2WsX171aJLFXny5MmbN29KShAZCUdHHK8cbZWairAwWFhg8OBqy6VqD/n5iIiArS3++ku1A2nz0mhUFDp2VN96mjxow9JoGUIhxo9H69YMBzFLRbmlUUNDZ6I7RGKi0VFRUapwTCq6pdHaIok+PHGCCgrI3Z2GDKFTpygtLS04eHSHDsO6dTvzww8UGMj8uF26kFhMCQnMW/4vMWQILVxIXbvSjRtSPrWzszMzM+NwKCqK1q6ljz+mqCgqO0Dq7V26WBoURNOnU7NmpZUNJAB49eoVtKMqy6lT5OlZOhHs21fT3miIw4dp9WraskXr1tO0Bz09WrGCBg2itm3pllbGFgqFYiIXIhaRb1raPU278x4qV2RVotIZYUUyM/Hzz3B0FBsaehOlEF0xNvbOUu6YaY388gv09HDhAuOGlUGbZ4QS/vgDjo6yK8C9eIFu3RAYWLXegoSKE8RNmx44O7cWCHo2bOjzUKO17QsKEBEBGxvsVKZ8ujJo54zw0SPY2NSUWkhTaNWMsIw//oBAAJXWKVduRujqGkIUShTDZtuffz/brMrQnRplnmfPXpqZ9ZRkO+NyJ546dYpZ+5cuQU8PCxYwa1V5tF8IAaxdC0dHyEx9LBYjJgY2Nti/X3qDrCwsWQJz84lECURgsRKGD5+koNe15eTJlKFDJ0RE/HD06Fs3NwweXLXwkErRQiEsLoa/v6pql9cS7RRCAAcPwtoa+/apyr4SQlhSggYN4OGxbsyYT69cuaIix6SiE0LmEYvFrq5tWKy/iRJtbb2ZTVVaUAAzM3TvzqDJ2lInhBDA77/DyQl378pumZgIe3uEh1ebkrRPn7FEZ4hAdMbYeOykSThwAO+ibFTLnTt3eDw/omQ2e5WBwaBt29QxaEW0UAgnTEDfvlq6iau1QgjgzBnY2Kgq+Y4SQjh4MExMNLPFq9sjZB4Wi5WY+OeYMUeGDYs7dGijKaMH2Dt2JAMD+vtvBk1+KIwdS7NmUefOdP++jJYdO9KFC3TvHgUE0D1p+xRRURMFgnGWllMFgnF//DHJ3p7mzSOBgIKDafly2faVo6CAkpMpMvJURsZQonZicZi5+YOhQ1UyVh1i2zZKSKDYWGKoFvoHhI8PJSXRggUUFaVpV4guXKCdO7V9i1eLXdNK7O3t1679lXGzUVF0/jxdu6bV94o28/nnVFBAgYGUmEjOzjW15HIpPp6WLCE/P/rtN6qiN56enjduHL527VqLFlHm5uZEFBFBb97Qv//SgQP0/fdkaUm9elGvXhQQQAYGyjt8+zadPk2nTtHp03T9OrVoQY0atTIyCs/NHcBiXXNw4Cpv+j9BWhpNmUL//ENmZpp2pW4iSUbfuzc9eUIrV2rywdKrF3XqRCEhGnNALtQwOVUdal4aVREpKWCzsWKFpv14j7qyNFrGokVo3BhPnqCoqEjm6k1qKho1QmioYvm6rl6tVPM9JqZS9XAAOTk56dLOsGdnIykJ0dEIDgaPBxsbBAcjMhIJCeVLrzt27Pb1DRk0KOzZs2cK+MQQ2rM0mpMDNzesX69pP2pEm5dGyyhLRs/g8r5CS6NTpqBePU2mx5JzaZQF7TgmrhybNm06cODA5s2bNe2I8uTnk7U1BQXRrl2aduU9CgsLORwORxMvkwDy8vJMFE+rs2gRRUV9qa+fxOEUx8REDRvWv4bG2dk0fjxdukTbt1OLFooNlJFBR47QoUMUH0+WlhQSQkFBdO3aurlzl7NY5r6+drt3r09LY507R+fOUXIy3bpFLVuStzd5e1OHDjKmrRpBKBQKhcJ69epp2hEaMoS4XFq+XNN+1IjSt6iaKS6mTz+lu3dp717iMrHQkJeXZ2RkxJJjwfruXWrShJYvp3HjGBhXOcRisUgk0tfXr7mZZoTw2bNnv/32W0ZGRvfu3Qe+lyuisLBw3759ycnJhYWFbdq0GTFiRHXP4v+AELZqRS9f0tOnxNa+7dq6KIRnzpwJDFycn7+ZKN/a2v/58wsyu2zYQNOn07ff0pQpyrgqEtGpU7R/P+3fj8uXPcXis0T6+vqf6et/6uzcrk0b8vMjPz/66CPS01PGvtrQEiFctIi2bqWkJDI01KwjMqgrQkhEAH33He3YQQcOkKNjba3JL4RNmpCxMV2Q/SeoQuQUQg08ffPz89u1a5edne3v7z99+vTVq1dXaXDkyJHFixfb2tp6eHgsWLBg9OjR6ndSPXz5Jd24QWfOaKMK1lFycnIMDGyIiMiopERPnjJMo0bR8eO0bh2NGkW5uQqPqKdH/v40dy6dOwcLCzaRHhEZGhrFxRVfu0br1lFYGLVsqe0qqCWcOkXz59P27dqugnULFouiomjcOPL3p4sX6e3bt8+fP1f1oDExdP8+7d+v6nEYQqXrs1JZt26dj4+P5P/x8fGurq5VVpwr7ktdvHhRT08vv5oV7jq9R3joENhsbN6saT+qp87tEQIoLCxs3jzAxCTcyGiAsXGE/KkJCgoQHo5GjV5/9FFPgcDbx6fn69evFR191qxoLrczjxfapk0vTf3olEbje4Tp6XBwqDbQU9uoE3uEVdi5E6ama8zNffj8rj16DBeJREoYkWePMCMD+vr49lulvGQU7Q2fSE5ODnyXmiwwMPDu3bvp6ekVG1Rci8vMzDQ2Njb8z70fvnlDISE0fDgNH65pV/5z6OmxWawHBgYZnTq17NmTTpyQq1e9erR4MTk6Rl+7NiY9PfX8+dHffDNf0aF//DEiNXVdQsL0lJS/NbKkXHcRi2nkSPr0U+rZU9Ou/Hfp3x96ekvfvEl++fLfM2dMk5OTVTRQz55kZ0c//KAi88yjNc0HIAAAIABJREFUgb/VFy9eNG3aVPJ/ExOT+vXrP3/+3Nra+v2WBQUFU6dOnTVrFrv6pcOzZ8+W7TLy+fxff2U+tkEV+PrW43JZq1YV5Odr2pXq0eweYX5+fg2/9+pITEx8/Lh1Ts6vRDhzpvUff/QbOLD+ypXF3bvLVapeT+85UTMiEoma79v39969RZ06iRRa1eTxeDwer6CGqiXaimSPUCxWa2UAAGvWbNi/P7m4uA3RuOnThdr8F1ERpW/RWiIUCiMjFxw9eiYw0Dcq6iuZu18VEYvF+vokWb0vKDDMzc3NV/zHLelSwx7hrl16qamGZ88W5udrvsqEWCyW50ekgWecoaFhcXGx5P8ASkpKpO7PFxcXDxw40N3d/auvvqrBmp2d3dB3sWAmJiYa3+qXh8mTWQ8fsu7dE2u/txoUQpFIpMTPx8LCgsV6SUREhSKRqEcPw7170bev4dy5GD1a9rmwr74ac/78uJycISYmcUOG/Dh/vuGYMTRgAEaOhL//fzyyWyOHZVavjp0zJzUn5xsWa83MmcuMjZU6sKQJlL5Fa8mPPy5avVqcn7/p1q0lZmar58yZrlD3CRNCly3rDNjk5+fu2xfUrZueone15LuuTggLC2ncOPbYsfDwqEWYLXNIVnFlt1Pl8qx0Jk+ePH78eMn/nz59ymKx3rx5U6VNcXFx3759Bw4cWPNGS53bI7x169Zff4HNVnlJHUaoi3uEYjGsrSeZmPiYm7c0M4s7cwYArl+Ho6O8KStv3bq1cePG2+8yeT98iJgYeHrC0RHh4Th3Tgmn6gYa2SMMCfmU6BIRiB61a9dfzaPXBk3tEXbtOoLoJhGIbtvYDFu9GlevQqHNvsePH1+7du31a3G7dhg1CnLsoFWi5j3Cnj1hYaGYPypFe/cI+/fvv3fv3pycHCLasmVL586dzczMiCglJeXGjRtEJBKJRo0aVVxcvHnz5v/MRsuFCxcMDOybNQvt18+jZ8/UD7akjqr5/XdydPztxYvjr1+f37ZtcEgIpaZS8+aUlERr1tDMmbItNG7ceOTIkY0aNZJ86ehIU6bQ+fO0fz9ZWNDgwfTRRxQVRXfvqvYb+UDo3NlXT28N0WMjo5Vduvhp2p06QL9+nTmceURnTUx+6tKlS2oqDRlCFhYUEEAzZ9LevZSZKcOCvb29m5ublRXr338pI4MGDCCmlqNPnKCDB2nnzjp4DF4NmlwFsVg8dOjQJk2a9OnTh8/nnz59WnK9R48es2fPBhAXF0dELVq08H7HoyrZO95Rh2aETZq0J4onAtFxGxtvTbsjF3VuRvj8OXi8SkWs9u2DQIDUVADIyICfH8aPr9XrqkiEpCSEh4PPh7c3YmLw4oXy1rQKjcwIV64U2dnF+Pn1mzNnft06Z6upGWFsrNjefsOgQV/8/vuGspnZmzdISEBkJIKCYGICNzeEhmLVKly9KiVleXFxcW5uruT/QiE++wy+vnj5Ul4HqpsRikSwskKfPsp9W6pCq6tPiMXic+fO7d+/PzMzs+zi48ePJQVRc3Jy7lamuu+kDgmhg0MbohNEILpiadlS0+7IRZ0TwsGDMWtW1Ys7d8LGBpcvA0BuLrp1Q79+qP0DXyhEQgJCQ2FuDn9/rFqF7Gzcv3+/T58x7dr1O3gwobYDqB31C+H9++Byod6yPIyhESHMyIC1NVJSampTUoKrV7FqFUJD4eZWmguwLJnfunVbeLyWPF7bwYPDJHomFiMyEs2bQ876m9UJ4ZgxMDautq6LpqiVEL569erPP/+cM2fOhAkTwsPDf/zxx4SEhDyFcjKqhTokhKNHHyBqyWKNY7MbrlmzVtPuyEXdEsJ9+9C4sXSF+/NP2NiUPnCLijB4MDp3BlNFtHJzsWULgoNhZgZT0wCiRKI7PF47zdb1VQI1C6FYjG7d8PPPahuQYTQihJ99hsmTFevy+DG2b8fUqWjTBkZGYn39VkQFRLC0HH2yQgHfJUtgZydXSXCpQnjlCthsbNmimG9qQEkh/Pfff/v06aP37sB4xRg+Y2Pjzz77TM1lFWumrghhcjLYbMya9Wjp0qW3bt3StDvyUoeEMDsbjo44dKjaBnFxsLXF1asAIBQiLAw+PgyXvX35Umhq2kZSt5nFmtOy5d9jx+LHH7FlC06dgrRE3NqFmoVwyRK0bQuhUG0DMoz6hTApCfb2ePtWeQu5uSJzcy8iMRE4nOlr1lRat9i5E9bWOHZMphEpQujggHbtlHdMdSgshPfv3+/atauenl63bt1Wrlx56dKlwsJCyUfZ2dknTpyYP3++p6ennp7eZ599VrbErFnqhBBmZMDICAMGaNoPxalDQjhpEsaOldFm+3bY2uLaNYn90uWganaflcTLq5ue3maiJB6v9b//voiLQ3Q0wsIQFAQXF9SrBxcXBAUhLAzR0YiLQ2pqeWJ+sVg8cuRkLreVnZ3n4cOJTLolH+oUwnv3wOcjLU09o6kENQthSQlatcKOHbW18+WXUVxub0vLMAeHLgJBUWgonj8v//TIEfB42L69JgvvC2FUFPT1kZFRW99UgcJCuG3btvHjx1d3LKWM1NTU4ODgNO24heuEEDo5oWFDTTuhFHVFCE+fhq0tKmw3V8v69bCxwY0bpV8uXgxnZ9y8qayX7/H69espU74dOnTCBWlrTG/e4MIF7NyJBQswYQJ69ECTJjA0hK0tAgLQteuhevU+JwLRS2fnNoz5JDdqE0KRCB06YPFiNQylQtQshD/9hB49mDF1/fr1kydPCoXCnBxERsLKChER5S9kly/D0RHLl1fbvYoQPnkCPT0sWMCMb4yj1YdlmEL7hbB3b9Svz/ASnNqoE0JYUgIPDxnvsBVZtw6Ojrhzp/TLjRshEMg4faBSRCI8eoTERHzxxXYOZx4RiMQCgZf6PVGbEP7yC/z9tSjUTDnUKYQPHoDHw927KjH+6BFCQ2Fvj1WrSleq799H06aIiJDevooQtmiBZs1U4hgjyCmE/5EoPe3k55/pn38oOZmZMmA6pDJ/PtnY0JAh8rYfM4YA6tyZjh4lFxcaOZLMzalPH9q0ibp1U6Wj1cBmk4MDOThQq1bdDhzo+vp1YUnJNUPDkKKi/2YFhrQ0io6mlJQ6GGqmOSZOpK++IhcXlRh3cKANG+jsWZo+nVatol9/pU6d6ORJCg6mMWNozZqaqtv//jvduEG3b6vEMbUiUypPnDgxYsSIoHcwIdKMoc0zwqQksNl1e/1H+2eEt27Bygr37ilsf80aODmVdzx2DHw+4uIUtsMsOTk5f/31V0rK6cGD0b8/1PyzV8OMUChEmzZYsUKlg6gJtc0It21DixYK539Rjvh4uLoiKAhXryI3Fz17ok8fVAkXKJsR5uTA0BBffaUOx5SGmaXRkpKS5s2b79mzJ/UdDLnHDForhOnpqFcPQ4dq2o/aoeVCKBajSxflXzVWr4azM+7fL/3yyhU4OGD1agAoOyamKYqK0LMnRo1S6/qhGoTwp5/QpYuUEO+6iHqE8O1b2NsjOVnV45RTXIyYGHC5CAvD06cYOxZt2lTa3CkTwg4dIBCozzHlYCbFWkFBQevWrfv06VOW5EU989Q6jVhMPj5kZ0fbtmnalf80a9dSdjZNnKhk988/p6lTqVs3evqUiKhFCzp2jObNe8bjtXd07Oz2f/bOO66p6w3jT0LYoOwNgigqiovhQEUQN+7iRqs/xS21teJqocMKTrBqxY3gQlssqKi4AFcVUXAPEEEFEUFAdpL390dQUZkhyWXk++Hj53Jz7jlP5CTPPee+5z2WfSSweWllyMnh6FEkJ2PhQqYkiJ779+Hnh927G3nuctGybBmGDkXPnpJrUVYWHh549Ajq6ujaFa1awdkZffogJeWzYidP4tIlhIdLTphYqcYIVVVV1dXVnz9/Lhk1jYPBg5GVhbg4pnU0al6/xvLlCAio087vHh6YMweOjnj1CgDMzGBtvSYzc2lGxuWHDz2XLVsjKrVCoKSE8HD89x9++olBFSKDy8XUqVi9GiYmTEtpONy4gdBQrF7NQNMaGvDxQXQ04uIQFARra/Tsifj4sle5XEyYgAkTYGvLgDZxUJURmpuba2ho7Nmzx9TUVOMDElPWQPH1xblzOH8ezZoxLaVRs3AhZsxAly51rWfRIsyeDUdHCIZ/RUW5gAEAIoOMjJw6y6wTzZohIgKhofCt9Q7B9Y5Vq6CujmnTmNbRcOByMWsW1q2DujpjGiwsEBKCvXtx/z4UFeHkhJkz1+jotG/WbBZQuG8fY8JETlVRo4nSBPu1JDoay5fDz6/x3CjVTyIicOsW9u4VTW3ffw8+H46OuHABy5e7//ffzKKiIaWlJ58+3ZKSwvAIRksL587BwQHNmmHOHCaV1IXbt7FlC+LipJOitcDfH82bY8IEpnUAjo6IjcXRo5g588DOnVsAP+BecbE1m32faWkioxYhzDweryHuuy0xXr/GwIGYOBELFjAtpVGTn49587B1KxQVRVbn4sUYPx5OTjAz656QcPzQoe7JySfmzetmb48bN0TWinDo6uLMGaxZg927GVYiHCUlmDoVGzfCyIhpKQ2HlBT4+GDHjvpy68BiwdUVxsa+wHxgFLCSSERbN9UPPjPC7t27b9myRXBMRJMmTbp27drHVw8fPqykpCRRdQ0HPh82NjAyQlAQ01IaO8uXw9ERzs4irtbbG66uGDgQaWnFJSUlRUVFHh7Yvh0uLti/X8Rt1RYTE5w5g59+wpEjDCsRgl9+gakpJk1iWkeDYsECeHjgw56Y9QUHB8FecvnATaCUaTmi5LOp0bS0tNzcXMExER04cGDYsGHdu0t3y6ye/v2Rk4N795jW0di5cQOHD+PuXbFU/uuvSE6OsbNbLCc3TEnpj7//3jB4cK+zZzF8OO7cwerVTN6et26NiAgMGAAVFQwezJiM2hIXh927cesW0zoaFKGhePgQISFM6/iKzZs3h4d3TklpzWLJ+Pp+x7QcUSLN7iAC/vgDUVG4cEEaICNeSksxYwY2bhRjpp63b3dzudsLClZmZu7w9d0FwMoKN27gyhWMHw9mnwx07Ihjx/Dtt4iOZlJGzSkuxpQp8PeHnh7TUhoOeXnw8MBff9XTvELPn99+//4Jj5fy448/Mq1FlEiNsK4I5qw2bYJ0jaW4WbcORkbiDR/Q1VVns5MBAM9kZMpipLW0EBkJeXk4OiI9XYytV0v37jhwAK6uiI1lUkYNWbECHTrUIvudFAA//wxnZzg5Ma2jiSE1wjqRlobhwzF5MubOZVpKY+fJE6xfj82bxduKj8/S9u39tbVtzcz+vHbNc+fOsvPy8ggMxKhRsLPDzZvi1VA1/fphxw4MG4b79Ttk7+pV7N+PTZuY1tGgSEjAwYONYbVMg0OadFt4BBlkLCwQGMi0lMYOEebOxcqVMDMTb0M6OjoJCeeJiMViPXmCkSNx4wb+/BNycmCx4OmJVq0wZAgCAjBypHiVVMHw4cjLw6BBuHhRXImY60hBAb79Fps3Q0eHaSkNBz4fs2bBxwfa2kxLaXp8OSJcsWKFjIyMjIyMrKwsgIkTJ8p8YPLkyUworL84OSE3F1euMK2jCbBnD3JyJLcuhcViAWjdGteu4c0bODl9mhEdMwanT+O77+DtLSExFTJpEn75BU5OXya+qicsXYpu3TBmDNM6GhRbtkBODlOnMq2jSfLZiHDs2LGZmZlMSWlYLF+OS5dw4wZUVJiW0tjJzMSKFTh5sk7Z1IRDVRV//401a2Bnh6NHYWcHAJ0748oVjByJx4+xezcUFCStSsC0acjNRf/+iI6Gri4zGirk0iWEhn5KxyWlJqSl4fffcfFifVk42NT4zAjXrl3LlI6Ggp/flqVL13G5Mnz+6l27XOue4ktKtSxYgG+/FUE2NeEQzIi2b4/hw7F6dVmSMAMDREfjf/9Dv34IDWVsAtDDA2/fYuBAXLjAZCKu8uTnY9o0BARAmo2xVixYgFmz0K4d0zqaKtJgmVpQUlLyww+ri4uv8nhxRL/Z2opnOZuUckRE4L//sHIlwzJcXBAdjbVrMWsWSksBQEEBwcHo3x89ejC5fvTXX9G/P4YMwfv3jGkozw8/wMEBQ4YwraNBceoU4uOxfDnTOmpAUVHRqlWrTp48ybQQEfPJCJ8+ffrmzZuaXPPgwYPs7GyxSaq//P13Cp9vAegBzYBOUVFRTCtq5BQUYP58BARAWZlpKYCFBa5dw+vXcHZGRgYAsFjw9sYff8DJicn9aNasgZUVRo5EVlbBw4cPS0sZS/lx7hxOnMC6dUy13yApKMC8efjrL8bm2GvOy5cvFRXNVq++7eKy2Ni4K9NyRMknI4yPj2/ZsuWPP/54r5L7WyKKioqaOHFi586dc3IYTswvSdLT4e+PTp2wcqUZkA1sBHYB59pJJzLEzPLl6NMH/fszreMDzZohNBRDhsDW9tNKvnHj8O+/mD2bsah3Fgt//QUZmZsGBr0cHH5v1apHamqqJAXk5eV5ePzk6Dhp/Ph/du6EmpokG2/w/PILevQQfcpAcTB27FhgEnASiH/xIoNpOaLk0zPCMWPGKCsrL1myZN26dW3btu3evbuFhYWGhgaXy83KyoqPj7969eqrV68GDRp08+ZNU1NT5jRLCB4PFy5g+3acPo2BA7FuHXbtkklJ+Y3HWwOQrKyGpaU0+ZxYOHfu3OjR8woKuDIyi1NTZzMt5zMEjwzNzDBoEDZswJQpANC9O65fx/DhuHUr6d07r+zsdz//PG/o0EESUyUjg+LiNcXFQRkZ7Vms8NWrt27dKrld7KZM+e7ECbvS0kny8t+rqOgC9hJruqFz9y5270ZCAtM6qiM/H6GhuHtXB+ACAPgMCxI1nwXLDBo0aODAgRcvXgwMDDx//vzeD/vcsFis9u3bjx07dsaMGe3bt2dApmR5+BB792LvXhgZwd0du3dDWRmzZuHcOQwc2LNDh4FFRUXFxTNHjVI6fVqaVk30DB06s7j4CNCWxxseGdls4sSJTCv6krFj0bYtRo3C5cvYsgUcDgwNERUFff0p79+vAYymTZt440a7Fi1aSEzSx2hDIkmngrt5M6G0dBeA4uKpFy5ctbeXGmGNECwcXLUK+vpMS6kEPh9XriAoCIcPo1kzsNm7WKz2RLeBHHn5mUyrEylUOe/evXv06FFycnJ+fn4VxRgkKCho4sSJoqotJ4cCA8nZmQwMaOFCSkj49NLixaSuTsOHU0lJ2Rk+n+bPpx49KDdXVO3XOwoLC0tLSyXfLpttBhBAwNo5c+ZIXkANycwkZ2caNIiysoiIuFyujk43gXJFxZ+PHz8uSTE3bsTq6XXR0ZmkpmZtYJB67VqtaygtLS0sLKztVbm5ZGAwncXaCtxXVx906dKlWjfcMOHz+Xl5ecJd++TJk++++2nkyA02Nnk8nmh1iYa7d8nTk/T1ydqaJk0iHR1yd6d374iINm7cGBp6k8OhpUuZVlkDeDxeycdv7cqpygjrP0Ib4dOnT3fv3n3jxg3Br7Gx5O5OGhrk4kIhIfTFl7+3NzVrRgMGUFHRZ+f5fJo7l+ztSdiPQ32HKSPU1u4IeAOhbLbp9evXJS+g5nC55OlJrVrRnTtERNbWAzmcvcBZDsd6+PDXGRkSFVNQUPDgwYOSkpITJ0hXl3burN3lQhjho0dkaUnTpuUuWPCTs/PEI0dCa9dkQ0ZoI8zMzNTX7wocY7H+7NjRReTC6kJqKvn5UefOZGJCnp70779kb0/W1lT+U/j+/Xs+n797N7HZFBXFnNaaITXCSrl5M05Ly0ZGZrOa2sAxY/a1akXt2pGPD71+XUHhTZtIRYX69fvSBQXw+TR7NvXq1Ti9kCkjXLkyx9h4cpcu/U6cOCH51oXg4EHS0aEjRygrK+vHH3/59ttFsbF3vLxIR4cCApiR9PgxWVqSuzvV4EugjNoaYXi4MHbbaBDaCE+ePKmsvFIwc6Cj06MmX9PiJju7bDJMU5Pc3SkmhvLyyNOTtLXJz4++GLMKjJCIXFxIWbm+f/VJjbBSZs9eCpwBCMjS1HT+MCysgMBAUlKivn0rdkEBfD65u5OzMxUU1FZIfYcRI8zKIh0devxY+HknRrh1i8zMaPbsZ8OHT+/Zc9SZM2eJ6PZtsrGhQYMoOZkBSbm5NGoU9epFaWk1Kl9zI+TzyceHTEzov//qpLBBI5wRXrtGvXsnysj0AbKABDOzbuLQVgVPnjzp339C584D/vknrKiIwsLIze3TZFhxMRFRWBiZmJCbG1U4pfHRCHk80tMjW1vJvoFaIjXCinn/nvr29Wex1gLEYp0fMmRKZSX/+YcUFalXL6r2y4HHIzc36t+/+pINC0aM8Mcfac6cOj2AYYo3b0hZuTdwAXiird3j+fPnRFRaSn5+pK1NPj7E5UpaksCxjI1r5Fg1NMKcHBoxgvr0ofR0EShsuNS2iyYkkKsrmZiQnx/t23e0bds+9vYj7927Jz6FFdK6tT1wDUhTUOitppbYrx/t2UM5OWWvPn5MAwZQ58505UqlNXw0QiJKSiIOh5YtE79uYZEaYQWEhZGpKY0eXejkNFFPz7pLlwGpqakVljxzhhQUqEePmnobl0uTJtHAgY3KCyVvhK9ekaYmpaY2SCPkcrna2hUHyzx9Sv36UdeuFBfHgLDwcNLRod27qylWEyOMjydz89rNuDZWat5F794lV1fS1SUfH4a/H7hcrrp6WReVl/feuzf840v5+SSYzPfzq+aOrbwREtGuXcRmU0yM+FTXCakRfkZiIg0ZQhYWFBlZfeHLl0lBgWxta9druVyaMIEGDapqHrVhIXkjnDWLPD2J6haSxyBdugyQkdkPxMjK2vTokR4f/+klPp8CA0lXlzw9GeghDx9Su3bVGFi1RnjoEGlq0t69opfXEKlJF71/n9zcyiyQ8UcnmZm0cCHJyg6RkdkBnNbTs379ISxCMEJwda3RKP8LIySiwYNJRaWePiysoRE2/lyjJSXw9YWdHWxtcedO9Rkcbt9Gv36wtER0dO2SHsnIICgIamoYPx7MZblqwDx7hn/+wZIlTOuoA2fO7J8//8G4cQevXdsxe7buwIGYMgWCxIUsFqZMwe3bSEyElRUuXpSosDZtPuWHe/261pfzeFi6FEuW4NSperdP0O3bt7OysphW8SXPnmHWLDg6omVLPHkCT08oKjImprQU/v5o2xYAkpIOLF36ZsaMyHPn9uno6CQlYehQLF6MHTsQEiLkNibHj0NFpWEkx6kUCXiy+Kh2RHjuHLVtSy4u9Px5jSp88ICUlKhjR+Fv30pKaMQIGj36yzUYDREJjwgnTKBVq8qOG+iI8AuysmjhwrLQu/LTTWFhZGxM7u6SXoT6McilwjUplY0IMzOpf3/q27fi0AkGycvLU1VtxWb3YLNb/PDDSgm3XlkXTU4md3fS0SEvr7KFd8xy9iy1b0/OznT3LhHRq1evpkxZOHjwlIsXL/v4kKYmeXmVxcjUkK9HhET09ClxOLRiheh0i4j6PjX69u3byp7P1ZwqjPDlS3Jzo1atKCKiprUlJ5OKClla1nUSo7iYhg2jb75p8F4oSSNMSCBd3U/G0DiMUMCDBzRwIHXuTNHRn05mZ5O7OxkaUqjEl96FhZGOTgUznBUa4a1b1LIleXoyEOlTLd9+uwzwBAjIZbNbHzhAcXFUq+Qfe/fuDQoKEq71r7toSgotXEiamuTpSdnZwtUqSp48IVdXat2aQkI+nbSycmKxwoA4GZluzs4vajhCKE+FRkhEO3YQm02XL9dBsRiov0bI5/PnzZunoaHRokULOzu7jK/uM/l8/siRI42MjABUnaiiQiMsKSE/v7I7spo/5EtPJzU1atmS3r+v8TupnOJicnGhsWMbthdK0giHDqU///z0a2MyQgGCxzAuLp8tpYiKojZtyNVV0oOtBw+obdsvHxl+bYT795etj2Sc4mI6f568vMjFhSwsSFWVWCxisfYBvwMEFAFWiookI0MAsdkkL0/q6mRqSl270pAh5O5Ov/5KISF069an0Y+mpiWLNYrFGq6r21EISeW76OvX5OlZZoFv34rqTQvP+/fk5UUaGuTl9dkD6exsrqrqx3iuX8LCwoSqvGIjJKIBA0hFpXb3IuKm1kb48OHDkPJ3DmLj1KlTxsbGmZmZRDR+/PgFCxZ8XWb37t0JCQnq6uq1NcKLF6l9e3JxoaSkWkjKziYtLTI1FeXz3qIiGjKEpk79cjlqA0JiRnjpErVo8dkntvEZIREVFJCPz5e3aAUF5OlJWloUEECpqalr1/odOnSYK/7x19u3NGAAOTh8yiNR3ghLS8nTk1q3LsuYIz5ev36tpdVBRsZETa1t0ocPbXY2hYSQhwc5OJCxMSkqEkAcDmlpkbU1ubnR5s30+DG9fv1aXt6YzXZlszuMGzfzY5337tGxY7R2Lc2ZQyNGkJ0dmZuTtjYpKBCbXeaUcnKPgaEfMvk5Ojo+/N//yNOTNmygAwcoKoqSkqr65Lq4jAH0ARMZGfMlS3iamrRwYU3Xa4oVPp9CQqhFC3J1pZSUT+cFybOaNycVFUcWKwy4raPTTbg5uSqMkMcjXV3q0UM47WKh1ka4e/duS0tLwbGuru5lsQ1xp0yZsmTJEsFxTEyMhoZGZSU1NDRqboRpaeTmRkZGFBhYOz3v35O+PhkYiD7qqaCAnJzo228bqhdKzAj79v3yr9YojVBAaiq5uZGx8Wdv+fZt6tQpTV6+M4u1W0Xl+3HjZktAiSA/nIkJxcYSlTPCN2/IyYmGDJHE/F63boOALQABBxQUJmhrk6wsAaSgQEZG1Ls3LVhABw5UOswqLCzcv3//rVu3at4ij0d37tDGjQ9YLGuAD/CBLpaWr8zNSV+f1NRIUZE4HGKxCCAWizgcUlSk5s1JT4/MzalTJ+rThwAD4A1AwKwWLXzqyZLK69epe3fdewk3AAAgAElEQVSys6OPmWbfvqWAAOrYkdq0IR8fSk8ve0Y4aJBbTIyQWWGrMEIievCAZGRopaSf2FZKrY0wNDRUT0+Px+ORmI2wT58+27dvFxynpaUByK0kZqAmRtiyZddDh0KWL09WUyudNCmztl+eRUVkYkI6OuIKW8jPJ0dHmj69QXqhZIzwxAlq0+bLOeRGbIQCzp8nKytycvo05AoKOiAjs14wRtHWtpGYEsGiiD/+uNOxo1Pr1r22br1oakqenpLosdnZpKQ0CLgJEJDE4fReuZLOnq1d7IbQ9Ow5nM02Y7NN+/b9prIyyckUE0OHDtHGjbRsGc2YQaNGUd++BBgBfICAVWz2vI4dydWVfvqJDhygmzdF83ilVqSn0/TppK9Pu3cTj0c8HkVGkqsrNWtGrq4UGUmVO1etqdoIiWjTJmKzy+6uGKeGRvhpGyYbG5t3794NGDDA3Nw8Nzd3zZo1uhXF0gYEBNQxTjUvL0/xQyixsrIygJycHFVVVeFqS05eOWHCL0pKa4yNv3/3rpjoYF5e9VcVFRXFxsa2a2fl5GTw/j0rPj4foJpcKASHDrHGjFGcPp2/aVPRx71yGgRFRUUcDofD4VRfVFiIsGyZkpdXSWEh9/PzVFBQQETia5pZbGwQFYWgIFlnZ/lRo0pXriwxMNBt1uxodvZCICkrS9HRkTdtWomLC1dWVrxKhgxBeDjb3n4I0Vqg+dy5323ZEubmppafL8ZG799nL16scOWKjILCTMAdmM9iBcyY0WPJkjwAxcUoLhZj6wJOnQrOy8tjsVgqKip5lXz4NTSgoYFOnb48r6cnV1AwAOgC7Dt16pCeXv6DB+yHD9lhYeyHD9mPH7Pl5NC2Lb9dO37btvy2bXlmZtSiBf/jx//48VOrV2/X0tLYsGGZubl5Xd5FaSl27pRbs0bO1bX0xo2SN29YK1fKBgdzNDVp2rRSf3+uigoBeP++Lo18RkFBAY/HY1X+Xfbtt/jnH0UnJ5mnT9/XagWaOODz+fLy8rLVforKu2J4eHjv3r319fXZbLaqqqp6RdTdogcNGvTnh7iI5ORkFotVXMkdYE1GhBxOkaLiyoiax4YS3blzh8MxZLMHAG2VlM5lZtb8UiF5/5569nyrozNQR8e6U6d+afXhYUINkMCI8OBBsrWt4Ha10Y8IP/L2LS1cWJbR49tvFygqttTS6nDhwuWQEHJ2Jj098vSkxETxakhKSmKzbQSDURZrxscJG3Gwfz+1a0dsNrVrR4LcOyEhIaNGjd23b5/4GhUH7u7uPXv2rDBHGo9HSUkUEUHr15O7O/XpQ9rapK5O3bvT9Om0dOmzZs3sgTTgv1ateta2XS6Xu3Klj53dsBUrVh8/XtqmDQ0eTPHxJOgwgqzZ5TM5iJxqR4REVFpKmppkby9GGTWkTlGjYp0a9fT0nDx5suD48OHD7du3r6xkTYxQRSVOW7vn89pEAXfrNgg4ABDwRFm5bc0vrAvu7stZrAMAsVjHJ02aL5lG64i4jZDLpbZt6ezZCl5qOkYoIC6ObGyecDg9gVg2e3efPqMF5x8/Jk9P0tEhZ2cKCRFXYjMejycjYwjcBBLZbPNHjx6JvInSUlq5kjQ0iMOhwYNrF8tWP6ltF83MpEuXaMcOGjPmJIfz04c4nR6qqiX6+tSqFVlbk6MjDRtG48eTuzv9+CN5e9O6dRQQQIcP0/HjdPEixcbS99+vl5X9HnjJYi3R1NywcWPZWlWx9pDy1MQIiejePZKR+bQymCnqZIRhYWGZVQ6Unjx58tdffwmn7OnTp6qqqkFBQVeuXLGwsPhYz5gxYw4ePCg4PnXqVEhIiIqKym+//RYSElJZhwsKCjI0tDlx4lTNW3/7ltTUpgJHBI8llJUthHsXtWX0aHfgOkDAo759x0mm0ToibiMMCCBn54pfarhGmJ6e/vjxYyEu3Lr1LxZrneD7kcPpePr0pzD0oqKy+319ffL0pGfPRKi3jLCwMHX19ioqbTZu3CTamlNSaNQo4nCoWTNatKjxJOMVuotmZGTo6VkDp2Rld9raDs7JoVev6PFjio2lc+coLIwOHKCAAPL1pZ9/pu+/J3d3cnWlIUOoTx/q2pVkZccBDwX38fLyvVq3Ji8viW5vUkMjJCJ/f2Kzmcmv+xHxriM8fvx4mzZthLuWiM6fPz906NDevXv7+fl9/D9dsWLF2Q+jg0WLFrmWo7K5xFol3c7Pp8mTSUaGmjV7xma3YLNHs9ktNm8W0s5rS0zMZS0tW3n5tRxOz4kTz0im0ToiViMsLCRjY7p6teJXG6gR/v67n7Z2L23t0Q4Oo75eAsHnU1YWJSXRrVsUFUVhYRQURJs306pVtGQJde++FrADEoFQNtvGyYmUlcnamjw9KTKybG3Jw4dlu8QJbv9F+8cRbof6Kjhzhjp1IhaLjIxo1y4RVlwvqEsXvX///owZi1euXP2u9rln2rTpBUwHEoAZhoY9RRgFU0NqboRE5OhIamoSin6qkBoaIYuEikc4ceLEDz/88PDhQ2EeX4qO4ODgiIiI/fv3V12My8WSJdi8Gaqq2LABU6fi/fv3Z8+e7datm76+vmSkAnj+/Pm1a9datrSZPt183DisXCmxloVErMEy69bh2jUcPVrxq0SUn5+voqIijqbFBJ/P19e3zsiIBWTk5ef16jVOXr5PTg4+/uTmQl0dzZujWTM0b/7lT2ZmzObN3kVFaoCavHyCgsKNPn1gbg4WC3fu4MYN2NrC2Rn29ujUCadOYft2PHiAyZMxdy5MTACgpKTk9evXBgYGMjIyQujncrlcLldBqPCGoqKi1NTU1q1bA+DzsW4dNmzAmzews8PmzbC2FqLK+g4jXTQ+HsOGbXv58h8+v7m8/LtVq0b98MNcSQoAkJ+fr6SkVEWwTHm4XOjpoX17REWJW1fF8Pl8Ho9XbbCMGAMC6wN8Pry8sH49OBysXo0ffig7r6KiMnLkSAmLadGiRYsWLQCcOwdHR3A4WLpUwhLqC3l52LABZ88yrUN0FBYiOBhZWYIHP5CR4dnbs+zsvnS7KunN5fY6dChMVVU5OPgvc3OcPo2ICJw+DT09TJsGIyO8fo3vvsOTJ7Czg7MzZs1CbCysrdG5MwYNurN27eSSEqNmzV5fvRouyTu8jRs3LV68FtBTUMh2db156FBzNhtjxuDPP6GmJjEVjZzCQvj64q+/8Ouvs4qL2RERlwcP7j9//gymdVUDh4PoaHTsCF9feHoyraYKhBtv1nFqVFRUPTXq60vKyqSgQIsWSVJUjUhPp7Ztac0apnVUifimRn/6iaZNq6pAA5oaTU8v28jN2ZmmTduqrd1dR2e4s/NYnogW4vF4FBtLPj5kb0+qquTiQuvW0V9/kbs7WVqSjg6NGkWTJ5OCgitwGyDg74kT5wrRkNBTo3JypsBrgIA1SkprtmwRoo6GhyS7aHR0WSq+j2mAGKRWU6MCNm4kNpsiIl7GxMSI6nNRQ2q9jrAxERCApUtRWIi5c7FuHdj1b7MpXV2cOYO+faGsjLmSnttgmMxMbN2KGzeY1lFnbt3Ctm345x+MHo2LF9GuHYA5WVnj8vLyBEN/kcBmw9oa1tbw9ERmJi5cwNmzCAuDhgaGDYOVFfLzERWF4uK7gGBWU/HQIYXsbPTtiz59YGMDccxtFxRg3z788w/i4lBSogcoAwCa9e9/san1Z7Hy7h08PXHiBP78E6NGMa1GWL77DuvXbx48eCebrSUv//TVq9tq9W2uQDibrbcjwt27SVeXOByaPJnJJ7Q15PlzMjOjbduY1lEJYhoReniQh0c1ZerziJDHo7AwcnYmAwPy8mImyTKX++UwUUvLBugAuANtVFS+U1AgAwPS0SEFBerShZYs+RRxUxnVjggLCykoiAYPJl1dYrFIXp7ataNFi2jcuB/YbEsWaxqHYyhcxGxDRAJdlKntuqpGiBEhEXE4JkAeQMAvc+bME4ewCmkqI0Iejyc4CA/H3LlIS8Pw4di7F82aMaurRpiY4MwZODlBSQlubkyrkQgpKQgOxr17TOsQipwc7NqFzZthYIDvvsOoURAqKkUEyMh8GiampSEiAjExKkApcANQ0dfPnzkTXC6yspCYiLg4+Ptj2zYUFaFFC/Tpg3Hj0Lv3l/tOl5SUZGRkmJqafn4SISE4cABxccjIgJwcWrbExImYNw/lMqKs++GHcZcuXZo2bUO9u9NvmKSnY/58JCRg3z707cu0GtFAgGBqjlNaymNYy1cIaYRsNlusmbdqzuHD90+csDM2vv7oEQYORHw8NDSY1lQbWrXC2bPo1w8yMpg4kWk14sfbG3PnCrkRNoMkJmLTJuzfj0GDEBICGxumBZVDXx/Tp+PCBauDB+15vGGysv6tWim+eIHkZKSmIiUFBQUwNUXz5iDC27cICcHevSCCpibs7DBuHL75Bjo6Ld+/LwFk2GwqLEypgfl9Ijr60rhx8/l89cOHIy9e/Ee4uFMpAogQFIQlS/Dtt9i/H/LyTAsSEdOnf7NzZ2ciS6KswYNDmZbzJbVYPiFIMSd0UlBxEBwc7OY2GRhrbj7x8uWRDe7r9SN372LAAGzejNGjmZZSDpEvn3j8GL164dEjqKtXU5KYWz7B5XLT09P19PQEb/zCBfj74+pVzJiBefNgYCB5RTUiKytr3Li5Dx8m2tvb7NvnLycn9/GlggIkJyMlpcwXnz/H06dITERmJths8HggygS6Ak8BWWAAsFVOrrWREdq1g709NDWradrLyyE9/TCgJy/vu3mzzowZ08T7VusH4uiiiYlwd0dhIXbsQPv2IqxYlNRq+UR5EhMT79+/f+DAkH/+kUlKgqGhONR9SZ2WT4wfP97Z2XnGjM9ic48ePbp48eLXr18L8V8gZrTy8k6vWzfSwQG9ejXIiO0OHXDyJAYPhqwshg1jWo3YWL4cS5ZU74IMkpiY6OQ0rrjYSE7uxYIFIQcOtCwpgYcHDh7Eh0Tx9RQNDY3IyEMVvqSkBEtLWFp+eZ7PR3o6kpMRGJiwfXszQOCdOpqaO52cfFksqKkhORnJydU0nZ9fAjQDUFys+euv+S9eYMgQ2NjUxyC1eguXi/XrsWYNlizB4sWMTbmLFXNzc3Nz82HD0LYt7OyQmlqfesjXjw1LSkrk5OROnjxJRM+fPxccENGrV68AJNWnLIFBQUEs1ncslt7Zsxk+PuTiQhoaZGlJ7u4UEkJv3jCtr5bcukV6emWZiOsDog2WuXGDDA1runs1U8EyY8fOAS4ABFzQ1597+rQo96+pzwAGwGhgBqD3tpbxP7t2BWtp9VJX9zAxsY6MzPDyImtr0tIiV1cKCKgX29WKAxF20Vu3yNqaBg+WaKY0oREuWKY8OTmkrEwuLqJSVBXCp1gTGN6dO3eIKDg4uHXr1oLzXC6XxWJdv35dtELrQlBQkKWlZfnPrSCazs+PXF1JXZ1atiR3dwoM/GyzZiLKzc2NiYlJryf7aZbj6lXS0qIP9x4MI1ojdHamgICaFmbECHk86tJlBnAVIODqyJEzJCyAWbp169a2bVvhtkZ5/vx5VFSUYOcsAenpFBhY9jEUJIqLiWlUdxV16aKZmZmbN287fDgkN7dUkDOv5h8Nxqm7ERLRtWvEZksiJbfwRpiZmQngypUrROTv76+pqSk4LzDIRHFvCVMbql5Qz+XS3bsUEECurqSpSS1bkpsbBQRQTEySgUFnDQ0PbW3biIh6l/bz8mXS0aELF5jWIVIjjIqi1q1rkRpf8kZ44QJ16UJWVgmamp21tNz19Drf+bhnbtNA5LlGP1RLMTHk6UnW1qSjQ66uFBhI2dmfCuTk5Fy5ciUrK0vkTYsVobtoTk6OsXFXDmeLgsJSZeWJrq6UkSFydWJEJEZIRH/9RWw2RUXVvaaqqFPSbQMDg9GjRyckJFhZWenq6h4+fJiIlixZoqmpWdnegYxQ86TbPB7dvk3+/jR6NCkq/gQcAwh40bXrYHGLFIKYGNLREXsXqRYRGmGPHvRhZ5EaIUkjfP6c3NzI2JgCA4nPp7y8vNjY2Hq7ilF8iMkIy5OYWHZj2rw5WVuTlxeFhNzT0+usqTlfV7frtWv/ibV10SJ0Fz116pSq6nLBNiNqat1r8jVdrxCVERLRiBGkqCjeZ1h1MsLdu3ez2WwA/fv337FjBwB1dXUAa9euFbXOOlGr3Sc+snz5KhmZQICAuzo6o8WxqU3dOXOGdHXpP0a/GURlhKGhZGVFtcqsJBkjzMsjLy/S1CQvr8azPZDQSMAIP5KfT+HhNGcOKSvPA6IAAhL69ZskmdZFgnBdlM+nX3+9x2YPAkqAV0ZGXcWhTayI0AiJyNyczMxEVVkF1GlB/bRp03r06JGSkuLo6CgrK6umpnbr1q2ePXsOHTpUUkE8YmTx4jmhocOzsg7IyGSMHRtoa4vRo/HLL9DTY1pZOfr3R2Aghg/H8eP1a9VabeHz4e0NH5/6FCEG8PkIDsbSpXB2xt279etP3xRQUoKLC1xcUFwsv2dPHhGAvKgo+eHD4eKCoUMlFFsvYR4/xty5yM62nDfPJTS0h7Ky0s6dfzItimGuX4exMcaNw+HDjOoQoxeLH+FGhALevHkjSP+amUmenqSpSZ6enz26qA+EhpK2Nt28SSK8Bas5IhkRBgZSr161vkqsI8Jz56hTJ3J0pFu3xNRCg0SSI8KPpKammpnZ6uoONTGxjo9/HBZG7u6kq0uWlmXxNZJN0VwLatVFCwrKMrP7+dFX+1Q2MEQ7IiSiy5eJzaZNIt4QugzxbsxbT6iLEX5BSgq5u5OODvn4ULnYN+b5+29SVd2irm6lrW21fr1E05LW3QhLSqhlS7p4sdYXiskIHz8mV1dq1YpCQkRed4OHESMkIj6fn5aWVn5TAkHst2AZhrY2ublRSAjl5EheWlXUvIueO0dt2pCLC6WmiluUJBC5ERKRry/JyFBsrGhrJaqxEdan6SpGMTZGQACionDzJtq0wfbt4HKZ1gQAcHB4CwRnZ9968yZu9eo92dnZTCuqBQEBaNcODg5M6wDevcPSpejRA5aWuHMHrq5MC5LyARaLpaenxy43dS7IpOrtjdhY/PcfevXCvn0wMECvXvD1xaNHn67l8Xh37tx5/fo1A7prQHo6pkyBuzv8/REeDiMjpgXVV5YsgZMT+vZFbi4zAqRG+Blt2yIkBEeP4vBhdOiAI0dQ4wx04iInJ0dRUR+QATglJXq5TPWU2kBEfn7bevceu2zZmp9/LmFWDJ+PffvQrh1evcK9e/D2/jLZtJT6jJkZ3N0RHo6MDHh6IikJTk4wN4eHB06cKOjUqV+/fmusrEZt376PaaWfwedj+3Z06AADA9y9i4EDmRZU7zl1Cmpq6NGDmdalRlgBdnY4dw6bN8PXF9264dw5JsWYmZm1bl2qprawWbP5PB5r8eIW794xqacm7Nmz/+efb1y65FNY+P7oUV8Jt56fnx8SEnL+/HkiOncOXbogMBCnTmHfvoaX7FvKR5SUMGwYAgKQmor9+6Gignnzwu7f7//mTdCbN+d//vnPejKFA+D2bfTsieBgREfDx0d641Uj2GzExiIxEVOnMtB6vdhBon7i7IwbN3D0KObOhYkJVq9mJnqTxWJFRx87f/48m83u3r3vsmXo0gXBwbC3Z0BMDYmMvJaX5w605PEWnDsn0X6dn5/fqZPTq1fDZGUjmzc/oqDw16pV0onQRgWbje7d0b07OnfmTJ1aXFgIgJuVxVJXR6dOsLGBrS1sbGBhAcknRS4owK+/Yvdu/P47Zs5kQECDRlcXYWEYPBh9+2KaZDO3VzMizM/P37lz59u3byWjpr7BYsHVFffuwc0Nrq7o3x9//hlhatrN2Nhm27ZAiclgs9nOzs5OTk5KSmx/f2zciDFj4O0NXr3b1asMMzN7Fmsr8FBJacOAAb0k2fTly5czMpwKC1fm5u7Izo6Lj+dKXbCxMmLEcCurWB2dYTo6DkFBy9LS4OMDMzOcPo3Ro6Gmhl694OGBfftw717Fzzh4PB6fzxeVnvBwWFqWzcC7u0tdUBgGDMCyZXB3x/37km246lia58+fA4iLixNJAI/IEWHUaLUUFNAff5Sw2R2Bd0Chllavly9fSqbpr0lNJQcHcnSkFy/E2IpwUaPHj5OWFnl47Orff9Lvv28QLu5UuKjRx49p5MhbbPZIgA9k6et3FqLpJgtTUaN1JDU19f3791+ff/eOYmLIz4/c3KhlS1JTI3t78vSksLCyPOC//bZRW7uTtnZHX98tQrRbvou+eEFjxpCFBZ09W4d30kAQR9ToF/TqRc2a1TRBf9WIZvmE1AjLk5mZqa09QJAbSUHBPSLihsSa/hoej/z8SF+fwsLE1YQQRnj4MOno0LVrdW26tkZ45w65uZGGBnl60oIFv+nodDI2tjl5st4lkq3PNFAjrCEvXtCxY7RiBQ0cSBoaZGSUIS9vD/AAbvPm3aKjs58/p5r3uOTk5K5d++vqWvv7b/XzI21t8vKioiJxvoF6gwSMsLSUdHSoqyiy7tQps4yUCtHU1GzVSjY/34vHU1VSujthQsexY+HtDX19BsSw2fDwQI8emDABR45g2zYoKTEgozyCXC2RkejYUXKN3r6NDRsQGYlZs5CUhObNAazctGml5BRIaQgYGsLQECNGlP168WLeiBG6xcVsAIWF2rNnv8/LU8vKQkkJNDSgoQF19U//fn1gazsoN/dHoKOHx/zOnXvFxFi1acPku2tkcDi4dg1t2mDmTOzYIZEWqxPE0dfXr3Z736bDxYv//PNPaGFh0dixZwoK5Navh5UVJk3C8uXMRCTa2eHWLcyZAzs7HDwIKysGNAj46y/4+uLCBbRuLaEWL12Cry/u3MGiRQgIqO8b50qpVzg4mHXuzL9zZy6LxevSRfHs2bIlfiUlyMpCVhaysz87ePToszO5ucXAdADAhG7dtrRps43B99IoMTNDaCiGD4ejIyZOFH97Ihh8MoeEp0YrJDWVFi4sy9DG4E4ygYGkpUV+fqKss+ZToz4+ZGZGItyhq+qp0ZgY6tePzMzIz6+pzEdJgMY9Nfo1PB7v4sWL0dHRQkz0KSu3AkKBJDa767Fjx8Qhr94iganRj/z4I3E49OiR8DVIM8tICCMj+PsjLg7Z2TA3x9KlyMlhQMaUKbh0CYGBGD0aWVkSbdrXF4GBiIlBy5bibYgI4eHo1g1z5mDKFDx+DA8PyMuLt1EpjRU2m+3g4NC7d29W7eM7r10LNTPzUVEZ6unpMuLjfKsUUbNmDTp3RpcuYcbGPRwcRgj2yhUHUiMUDSYmCAgos0MLC3h7M5ArqE0bXL0KCwt07YroaEm0SITvv0dICKKjxbtdAJ+P8HDY2ODnnzFvHuLjMWUKONIH3FIYokOHDomJV9PSrv/xxy9Ma2nkfP/9kYIC3xcv/KKje3fqJK4MPVIjFCWmpggIwLVrSEuDhQV8fVFYKFEB8vLw8cGOHZg0CUuXorRUjG3xeJg5E5cvIzISWlqirJnP5//888/Tp7vfvXu3tBT79sHSEr6+8PZGXBymTKlfOzpJkSJFfISEHAXcgW7A4vR0cU12Sb9RRI+ZGQICcP487t1Dy5bw9UVRkUQF9O+PmzeRkICePd926zbG0NBm6FC3goICETbB42H6dCQm4uxZaGiIsGIAaN3aftWqtCNHenXqNLJFi9fbt2PLFly6hGHDpIuUpUhpWgwa5MRi7QFeAUf4/K6OjhBHinWpEYoLS0vs24fISNy8CQsLrF6dbWs7RFfXpkuX/hJIlq+jg+PHwef/dv36+FevYs+e7bZqlb+oKi8pgasr3r5FRARUVUVVaxmlpUhOTifaCczl82f16eN/6RL69RNxK1KkSGkQzJo1a/r0jsrKjmZmG44cWfP6NQwMMGQI3r8XZSs1MsIXL16cPHny5cuXomy5adChQ9l2FgEB62/enJSRERsf77F48SoJNM1mQ1ExFbADUFLS7fjxFJFkyisowLBhkJNDaKgoswlnZWH/fowfD11dEPGBZIDPYkX17WsisjakSJHSANm5c9P794+Skq5+8435/fs4dgx37kBTE99/D1ElyKvYCMeNG+ft7S04vnjxooWFxdChQ83NzY8dOyaaZpsYdnbo2vUNURsARG0iI1+fPCmJ/Q7nzBmvoeEO7GnefJGJybjWrTF2LM6eFX5vqZwcDBgAHR0EB0Mki0uTkuDvj/79YWaGvXvRowfu3MG+favl5PrKyJja2WH27NkiaEaKFCmNhWHDkJqKVasQEIBmzbBpkygq/XpFRWlpqby8/OnTpwW/du/evWPHjjExMdOnTzcxMeFyucKv6RA19WEdYQ2Jjr6kpWUnL79eU9P+hx9OOzuTpia5uVFkJIl1Tc7169f9/TffuXOHiN69o4AAsrKiNm3Ix4cyM6u59ot1hFlZ1K0bzZ1L5fYSFwYul2JiyNOTLC1JR6ds//EvFg2KaYd6KVXQ1NYR1pEm20UluY6wang8WrSIZGVJQ4OOHKmsjLC5RgVPsJ4+fSo4ZrPZ+/fvJ6JXr14BSEpKqpN2kdKAjJCInj17tn///idPngh+TUkhPz/q0oWMjWnhQrp5U3JKYmPJ3Z3U1cnVlSIjKy1W3gjT08nKijw9a9pEfn5+bGzsu3fvPp7JyqKQEHJzI3V1srQkT0+Kian0JqDJfsswiNQIa0WT7aL1xwgF5OfT5MnEZpOREV2+/OWrwi+oFyRUKy4uBhAREUFE/fr1A6ChoQGgyW7JVHdMTU0nTpzYqlUrwa/GxvDwQFwcIiKgrg5XV7RvD29vJCWJXYm1NQICkJgIe3vMno327eHri+zsSsunpKB3b4wZAx+fGtWflJRkYdFr8OAdbdo4HT16dft2DBsGU1Ns3w5rayQk4N49+PigVy9pFKgUKVLqhJISgoLw6hWsrNC7N3r0QLgln5EAACAASURBVGpqrSupwAjV1dUNDQ137tyZk5Ozc+fOLl266OrqAhDsRKGtrV1n5Th27NjUqVMXLFjw8OHDCgs8ePBgwYIFU6dODQsLq3tz9RyB/z15goAAZGejZ0/06gV/f2RkfCqTm5ubl5cn2nbV1eHhgceP4e+PmzfRqhVmzcKtW58KvHv3jsvlPnsGR0fMnQsvr2oqLC3Fixe4ehWzZm17+fKPN2+2vX59YPLkDTdvYtYspKcjMhIeHjAyEu37kCJFSlNHVxcnT+L6deTkwNQUQ4Z8SmlSUlJS/fUVjhP37NnDZrMBsNnsv//+W3By/fr12tradR8UHz58WFdXNzg42NvbW1NTMz09/YsCaWlpGhoav/zyS3BwsK6u7pHKZn8b2tRoDeFyKTKS3NyoeXNydqbAQJo3z0tbu7uWVrcff/xNfO2mpZWlDLW2poAAUlJqBbQATDicLbt2fSqWk0P37lFkJO3dS7//TvPm0fDhZGND+vokJ0eGhtSjB5mbL2ex/gYIuNWv3wQhxDTZeScGkU6N1oom20Xr29To15w4QUZGxOFQhw5/sdmG69efr/aSSpNuJyQkBAYGxsfHfzxz8ODBgwcP1l2ljY3N7t27BcfDhw9ftWrVFwV+++23ESNGCI537dplZ2dXWVWN0gg/kptLgYHk6JjGYvUVbIKoptY7ISFDrI1yuRQeThYWAUBPgA8UA61cXXOcnKhdO1JWJlVVateOnJ1pyhRavpz+/JP+/ZeuX6eXLz8F0aSlpbVsaaerO9zQsIsgTqe2NNlvGQaRGmGtaLJdtP4boQB/fwKsgHc//lh9dF+l6RqtrKysPt/UZ/z48XUYuZZRWloaFxfn4OAg+LVv374XLlz4osy1a9ecnZ0Fxw4ODjNnzuRyuZyml1lSVRVTpqBv3xJraxVBstn8fOV+/Ury8tCyZQU/Fe5DdPz4yYsXbwwd2tfR0eHrV0tKkJKCZ8++/MnOjgOMARYgByh07PiqR49mBgYwNoaKSvXK9fT0njy5mpaWpqur2wT/cFKkSGGchQuxaFEun9/MwKD6xYaVfkm9ePFi586d9+7dKywsPH78OIDjx4+rqqp+9DDhyMjI4PP5mpqagl+1tbXT0tK+KJOenv6xgJaWFp/Pf/36tWElSZ1v3LgxevRowbGOjs66devqIq8eoqGh0a2b2n//jQX49vZ6wcHNi4vfp6Wxnz1jJSeznz1jRUWxk5PZjx+zZGVhZkampnxTU77g4Pr1fZs2nc3NnbJr16oNG962ajUoOZktuErwb3o6W02N9PXLrrK3p0mT+GZm1KzZyhYtbAB34K2MzLvvvjMCyhI51DyhQ/PmzYuETS5HRKLNCSelWrgfYFpIw6DJdtGCggI+ny/Elh2Sp2PH1vHxA0pLNwIdqi5ZsRHeuHFjwIABLBbL2Nj4Y5jo3bt3AwMDHzx4UBdlioqKKPf0sqioSOmrjdUVFRU/FhAEr35d5iOGhoYfh6oqKipVlGy4hIUF3r17l8VitW/fHoCSEtTVYWn5WRkievkSSUl49kwmKUnm2jUcPMiKjT1VWroJMH73znDu3ICuXUebmpKZGeztMXkyTE3JyIj/YcDG/jxyyuDt2/uTJk3S0NAICnousXda/u0QUaP8a9ZbBC6oIMJ0QY2aJttFBe+6QRjhzZunDx8+bGtbfZeu2Ahnz57dpUuX0NDQW7duTZ48WXDSxcVl2bJlGRkZOjo6QitTV1dXUlJKSUkRRJ+mpKR8PdQzNDRMSUkRHKekpCgrK6upqVVWoYGBwdixY4XW01Do2LFjtWVMTGBigr59P52ZM6fNrl0RpaXuCgonV6xos3IlgPLdt6qurKGhERoayuFw2Ezs9UBEbDabkaabLOwPMC2kYdBku6jgXTcIIwQwbtw4Ho9XbbEK/orZ2dlxcXG//fZb8+bNy7/bFi1aAKhjxlEWizV69OigoCAAhYWFR44cGTNmDICCgoKDBw8WFhYCGDNmzJEjRwTHQUFBY8aMaSj/6fWNNWtWDB8ea2JiO3782yVL5jMtR4oUKVLqIxWMCAWzkapfbSuQnZ0NoO6xD15eXk5OTgkJCa9evbKwsBDs7/z27duJEye+ePHC0NBw1KhRgYGBXbt21dfXT0xMPH/+fB1bbLKoqqoePbqdaRVSpEiRUq+pwNV0dHS0tbVPnjzZsWPH8kOxQ4cOKSsrW1hY1LHJVq1aPXr0KDY2VlVVtXPnzoKTBgYGycnJenp6ADgcTnh4eHx8fF5enq2trfShhRQpUqRIER8VGCGbzV60aJG3tzePxzM0NOTz+ffv3z98+LCPj8+iRYvk5eXr3qqiomLv3r3Ln5GRkRFMvQpgsVgfPVKKFClSpEgRHxXPc3p6emZmZnp7ewtiqdu3b89isaZOnfrrr79KVp4UKVKkSJEiXioOeWKz2evXr09MTNy7d6+Pj8+2bdvu3bu3Z88eOTk5CeurlmPHrm/bto1pFVKkNGwCAgLk5U2UlMznzp3LtBYpUiQNi4TepLUeEBwc7OamDbhfuxbSrVs3puU0NoqKijgcDiOpYYgoPz9fpSZpbKR84OTJiJiYmyNGOHfv3r2217JY+sB6QA2YeebM3v79+4tDYWOiyXbR/Pz8hrKOEACfz+fxeLLVbSNe8Xfc8+fPK1t70bJly7pKEzEDgYk+PntDQ6VGKKXp8uefO3766WxOzoTt25cfPLhiwIB++fl4/x55ecjJQW4u8vLKfn33Dnl5n37NzkZy8nnAEJgIABg5Zkzi//7X38kJ/fqh6a0Xl9IUqdgIu3XrJtie92vq3wgyFTh77NgRQ0O4u2PZMtS/6VspUsRIVhbi4rBmzb85ObsA3awsvWHD9nO5/RQUoKICVVWoqUFVFaqqZb+qq0NFBdraUFWFggLi4pCU1Bt4BSQAzYELenpzjxzB1q0oKYGCAnR00KYNrK0xYAB690aFEwRZWVmJiYnW1tZNcIG5lEZAxUa4Y8eO8lki8/LyoqKi/v333z/++ENSwmoKm+08YkQPdXXTw4exdi1++w0dOuDnn/Eh/6gUKY2NV68QF4dbt3DrFuLikJ2Nzp3RrFlrGZkzPJ6bgsKppUstVq6EjExVlWRkYM8ebN6Mli3BZsu2aNHj+fOhANTU7A0M2h8/DhUVFBQgMhJRUYiLQ3AwNm5EcTHk5aGtDXNzdO0KBwcMG4Y//vD18toKmMrJJT9/fqMumaekSGGEWjwj9PHxOXXq1MWLF8Wpp3YEBwdHRETs378fQFwc5s5FQQFYLNy9CwUFjB6NNWugr8+0ygaL9BmhhElKStq7N8TU1MDNbUL5pxqvXuHmzU8/hYVo3x7W1mU/7dqBzUZOTs7kyR7x8XednfsEBPhW8VDk8WNs2YKgILi44Jtv4OGBqVPh5QUej8vlcmVlFWbOxNOnOHECXyXVQFoaTp1CVBTu3MHz53j3Dnw+iOyACEATWD9+fMLBg4Fi+v+pVzTNLopG+oyw0v0Iv0aQ/zMxMbHml4ibL/Yj5PEoMJD09Gj+fPr1VzIyIhaLjIzI1/fTVnlSak5hYWFpaSkjTTfBzd5evnypq9uFxTqkpLRy4MD/hYWRlxe5uJC2Nunrk4sLeXpSYCDdvUvCbQbH51NkJLm4kJ4eeXnRmzd07hzp6NC+fWUFPu5HyOfT/Plka0tv31Zf7dOnJCPTGsgFCNjK4fgMGUJnzwqjsGHRBLuogIayH6EAHo9XUlJSbbFaTOjn5OQAqM87j7DZmDIF9+6hpAQBAVi1ComJ6NsX3t6Qk0OPHrh0qaxkQUHBiRMnMjIyGNUrRconIiIuvn3rRjSuoOC3s2fjg4OhqIgFC3D/Pl69Qng4fHwwZQrat0dt78WLirBvHzp0wJIlGDYMz57B2xsnT2LiRBw+DDe3L8uzWNi0Cfb26N8fH/aeqRRzcyxY8A2bbc1mu8rK/uHlNfnlSwwcCAUF9OiBrVvBr34zOClSmKZCe0xOTk4sx8OHD0+cOGFjY6OmplZcXCxqzxaeKnaoj44mKytycqIHDwQlqVMnYrNJQ4PGjHnI4Rix2cPYbJNt27ZLVHGDQjoilAClpXT8OI0fT6qq1+XkhgPFQEK7dg4iqTwtjby8SFubXFwoMvLTeT8/atmS7t//QsmXO9SvWEGdO1NGRvUNPXjwICQkJD8//0NVFBREDg4kJ0ccDnXqRFu2EENdSVw0nS76BY1yRFixEerq6n5tmQYGBuHh4aLWWSeqMEIiKi0lPz/S1CRPTxJ8wN++pTlziMNZBIQABCSqqLSRnNyGhtQIxcrdu+TpSXp6ZG1Nfn6UkUFr1mxu0cLO1nbIA8HtWy05ffrMTz+t/u+//4jo5k1ycyMNDXJ3p4cPP5UpLSV3d7K1pfT0Ly//2giJyMuL2rWjly+FkFPGsWPk4EBKSsRmU8uWtHIlNY4/bFPoohXSKI2w4mCZ8PDw8lGjHA7HyMjIysqqvuW/Lh8sUxlpafD0xJUr2LQJQ4YAQLdug65fnwaMAxJlZceVlMRKSG5DQxosIw5evsTRo9i7F1lZmDAB//sfWrcWQbWbN+/66adT796NU1HZbGHxc3a206xZcHeHuvqnMnl5GDsWCgrYv7+CBYKVbczr64s9e3DuHL7aObR2REfDzw/nziEvD2ZmGDECS5dCS4v/zTfToqPjbWwswsKC62HuqspoxF20app6sEw9pOoRYXnOnydLS3JxoeRkio+P53AM2WwXFqslh3PRyIgePRK30gaJdEQoQgoKKCSEXFxIU5Pc3CgyUsiYl8qwsRkKpAMEXOndewGX+2WBFy+oUydauJC+fklAhSNCAWvXkpkZJSWJRuqDB+TuTrq6xGKRvPwa4FvgATDfyWmMaBqQCI2vi9aQRjkibCqrXx0dcfs2nJ1hZ4e//+4wb94MPb2k0aP7JiV1MTJCu3ZYupRpiVIaIzwezp7FlCkwMsK+fXB1RWoq9u2Ds3OtY16qRlPTHDgLQFHx7NChrb5YRJiQAHt7TJ8Of/9q1hdWyOLFWLIETk5ITBSB1LZtERCA9HQ8fAg2OxJYBLQFlkRHp584IYL6pUipLZ9mvSIiInx8fKq9ICoqSpx6xIisLDw8MHw4Ro4Muns3i8+P/fffA0pKP1296r9jB+bPx/79OH9eNPNUUpom4eGn5sxZUVLCdXMbM336z0FBCAyEoSHc3LB+PbS1xdXu8eO4efMXW9sF6ekbHRx6fv/9Z7d1p0/DzQ2bN2PsWOGbmD0bMjLo0weRkbC0rKtgARYWcHExOnLEH1gJbFZSsh09GlwuTE0xYgQWLYKxsWgakiKlaj6NCGVlZVVqAINaRYKZGWxt7/L5owBFLnfszZsJAGbORFoaBEPD5cuZliilwTJr1rKXL8+/eXPb3z9+0KDbCgqIiUFsLDw8xOiC+/fjf/9DWJja9etBKSmxQUGbyj8R2bED336L0NA6uaCAmTOxdi0GDMDdu3Wt6iOHDu10cSls3nx4v37Jr1+vLi7G+fPo2RMHDsDEBM2bY8gQBAdL12BIETMSmKUVHzV/RlieM2fOqqv3B86xWLOnTFlT/qVt20hWloyM6OlT0alssEifEdaK1FSukpItQACpqHieOnVaAo1u3kyGhpSQUMFLfD55eVHr1vT4cY2qquIZYXkOHiQDA4qPr6XQ2vP2Lfn6UvfuJCdXFnG6aBGlpIi93RrSELuoSJA+I2wk9O/fLyTEc8aMU0uX9j5z5odHjz69NGsW0tNhZIQ2bfDTT8xJlNKgeP0aS5fCykrGwMC2efP/KSv/oqcX1bt3L3G36+uLzZtx5QqsrL58qbgYEyfi7FlcuSLi2f7x4+Hvj/79cf26KKv9Gg0NLFmCq1dRxTDR1nYwh2OiqGj+77//ileNlEZNVZHxGRkZSUlJ79+/L3/S2dlZzJIkgbNzP2fnfgBat8aQIbh6FR8TBWto4OpVbN2K777DgQO4dEmarVRKpaSkYP16HDiA/7N353Ex538cwF8zXbov3SWVM0QqV8iRI5Izt1wry6LF0jr2JyySY93kJvdNEnJsJeRYV4kVOdJ9KLpn5v37Y6y1KV0z851pvs+HP8b0PV7lq/d8P9/PMXIknj6FicnmP//8Mysrq3fvOWriXMGICLNmITISERFlNLpmZmLAAJia4soViGPE05AhUFWFhwfOnEHVlz6sDhcXuLgAQEYGNm7EqVOYMAFeXpeJdIDXfP6LYcPcCgv7SyIKq1Yq8z4xNTW1e/fuld+eKdVrGi1l/nxq04b+mRPjX4mJ1LQpKSrSqlU1PIOsYptGvyMhgWbMoLp1acYMSk6W6KmLimj4cHJxoZycf98sLi5etCigW7fh//vf9kaNaMaMKg/PqGTT6BehoVS3Ll27VrWziFDPntOBDcK2aMBeT48cHOiHH+jQIUmM2Zf+S1RM5Khp1NvbOzY29uDBg+7u7hMmTLh48eKMGTN0dHSCgoLEWZSZ8fvvaNgQ48aVfiBvZoanT7FkCebNg60tylmfkSV3EhIweTJat4aqKp4/x/r1MDaW3Nnz8zFgAPLzERoKLa1/31+wwD8gIPfatSVLl97s0OHw+vUiHp7xrd69ceIERozAlSviPVF5/P3Hc7mrgJ0czkQrK62lS1GvHq5dww8/QEsLysowNkaXLvjtN4SHl93d5syZM8HBwRIPzpI+ZZZQVVXVgwcPEtG4cePmz58vfH/dunXNmjWTqs8CIrkjJKKCAnJ2pgULyv7ql1vDtWuJiHJyctLT02t+UunH3hGWEh9P3t5kaEiLFlFWFgMBsrLI2ZnGjClj3s7Wrd2ANICA28OG/VSNg1f1jlAoMpL09ePq1XMxNm49cuRUvmTXeblz586AAUMXLlxY6rx5eXTmDM2cSS4uZGREXC5xOKSlRU2b0sCBtHkzpaeTsXFLLrc3l9vD3NyhGqeWzktUAmrlHWEZhTAlJQXAs2fPiGjy5MnTpk0Tvp+ZmQng70p2QZMIURVCIsrIoIYNaevWcjdYupQUFEhD439cbn0ut2GLFt1Ecl5pxhbCL2JiaMwYMjKiRYvowwdmMiQnf54apszfQj//vJTDWQA819Ly2rPnYDWOX71CSESNGvUCYgDS0Jh95MjRahxBAh49ouXLqW9fsrKiOnUISAD6CJtVOZzOffs++/lnWr6cdu2i0FB68uQ/zc7f8vQcBZgCFsrKlhL6BqRGrSyEZXSW0dHR4XA4wkWXzMzMIiIihO8XFBRAupdhqgl9fYSGolMn2NigR48yNli4EO7u2fb2p4B4QCE2tvPVq1fLe5LKqjWePMGqVQgLw+TJ+Pvv/7RGSlJCAnr1wsiR8PMre4PkZN82bdZraPgNGtR93LiRksxWVJQNNAGQl9cqIeG9JE9deXZ2sLP796+3buU4O6cLJ1omKnj2TOPRIxQUoLAQJSXg8T43pXK5UFSEsjLq1IGaGjQ0oKUFbW1cunQduA8YFxeP9PDwOHfuHDPfFUtEyiiEKioqzZo1u3v3bps2bdzc3BYvXhwQENC6des1a9Zoa2s3rL0zr9jY4OhRDB6Mq1fL6I8OwMAgn8vlCAQKAAQCs7i4XLYO1jL5+fnjx8+6detehw5Os2atXb1aNSICM2di27YyZqmWmJgY9OmD337DpEllb7BlC+LilG7f/kVVVbLJAABjxw7ZuHFkTk4nJaVdw4adYiBB1bVv37JlS93Hj20BgaNjw+joMiYUz8rC27d4/x5JSUhNRVoaMjLw4QOyswEoAMIleuzDwlR++AFt26JtWzRrVp0Z7FjMK/M+8cCBA+vWrRO+njFjhnCicXV19SNHjojyrrXGRNg0+sWRI2RlVcYiNUKWlm243N5c7hAFhY5cLn/gwNq2ytrX5LBp9JdfFisrrwOIy12npbVk0yaqVmOhKN2+TSYmdOxYuRs8fkwGBv9Za6l6qt00SkTh4eE7d+5u3jzx8OGaxpCkt2/fJiYmVmPHOnXqAx7AEsCsVatn27bRjBnk4ECamuTsTDNm0L59FBMj8rxSoVY2jVZqOER6evqdO3c+MPVspHziKIREtGgROTjQp09lf/X48eN79+7l8/nHjpGmJmlp0enTIo8gFeSwELZpMxqIAwiI6959lOQDlHL+PBkY0KXy56j5+JGaNKGgIBGcqyaFUOjOHTIyIvnoSUZDhw5t2bJlZGT0jBnUqBE9eUJElJtLkZG0bh2NGUP165OJCbm706JFdO4cZWb+u++zZ89+/NHXzy8g5/uPIqWSHBXC6n1KkjwxFUKBgLy8yN293AVrvuDzydubuFxq2ZKSkkQehGFyVQgfPSJPT9LXP66m5g6c0dXte/jwCUkG+NaBA2RkRDdvfm+bMWPI21s0p6t5ISSiGTNo/HiRxJF2X1+iwn+pb++G37+nc+do0SJydSV1dbK2pjFjaNmydD291sAFJaVtbdr0kXTuGpOjQmhkZOTg4BAYGJibmyvqYKIkpkJIRMXF1L07/fJLpTaOj6emTUlBgWbOFEcWxshJIXz8mDw9ycKC1q2jggK6evXa3LlLrl5lbqA4EX13EtEvdu6kZs3KmAuiekRSCD99ImtrunxZJImkWqlL9OlTatqUvL2pqKjs7UtK6MED2raNeva8wOX+JuywamDQvjK/pqWKHBXCwMBAR0dH4XNBLy+va9euSXh4UCWJrxASUWYmNW5MmzZVdvsNG0hFhYyMKDJSTIkkrdYXwidPyNOTjIzI35/5Z4FElJycvHnztnPnzq1YIWjShN68+d7GMTFkYECxsSI7u0gKIRFduEBWVuU+Wag1vr1Ec3NpyBBydKSEhO/t+OrVq7p1OwFpwH0FhXZbtlTc8iRV5KgQCsXFxS1atMjS0hKAubm5r6/vixcvRJdQBMRaCIno1SsyMaFz5yq7fV4eDRxIXC65uEhikidxq8WFUDgu0NiY/P0pP19856mC1NRUU9NWCgpblZV/0tefkZr6vY0LCqhlS9q7V5QBRFUIiWj4cPL1FcmRpFeZl6hAQOvWkYnJ9x7rEtGJE2ebN+/WpYtnaOizXr2oVStZ+vQsd4VQiMfjXbhwYfjw4XXq1OFwOKLIJjLiLoT0z/P/hw+rsMuff5KBASkry/wkpbWyEMbGfh4aLz0lUOjw4cMqKquELWZGRhXMdTJhAg0fLuIAIiyE6elkYkL374vkYFLqO5doeDiZm9OiRVTJdrRz58jSktzdK2gDkBK1shBWvAyTgoKChYWFmZmZlpYWCQegyhMnJ2zcCHd3JCZWdhcXF6SlYe5czJsHGxs8eybOfKxKi4uDlxe6dUOzZnj1Cr6+YGTUXXmMjMyJ7gMCIEFL63trRhw9ishIbN8usWhVVrculi/HhAng8ZiOwoTOnXH/PiIj4eEhHHRYgX798PQpHBzQujX8/FBUJP6IrP/6XiHMyMjYsGGDg4NDixYtduzY0a9fv8jISIklkx6enpg6FR4e+O+CVBVYuhTJybCwQLNmGDQIPB54PN7Dhw+Li4vFlpRVtlevMHkyunaFtTX+/hu+vkyOji9TcTHWr+9Yr15zExOnpk3HHzmysbwt4+MxfTqOHYOmpiQDVtm4cTA0xIYNTOdgiKEhLl1C8+awt6/Uwo1qavDzQ3Q0nj5Fixa4cEH8EVlfK/M+MTg4eMCAAcrKylwut3v37kFBQXmi6pomUhJoGv3ixx/Jza06w+ePHydNTVJXj1dQsORyOyoqml2WkU51Mto0mpWVtXv3njNnzvD5/Fev/p0mW/rGwX6Wl0c9e9KwYVRhE05hIdnb07ZtYokhwqZRoYQEMjCg+HgRHlKKVPISPXOGjIwoMLAKRw4Lo6ZNyd2dXr2qfjzxqZVNo+UOn6hXr56vr+/Lly9FHeyz5OTk58+f1/AHKslCWFxMPXqQtzclJSW9fv26SvuWlJCm5kDgAkBAlLFxazGFFC1ZLIS5ubmWlg5KSn+oq8+0tBwn5SWQiD59om7daOTISn3GmjKFPD3FlUTkhZCI/P2pW7cqr4woEyp/iT5/Ti1a0JgxVRjoUlxM69aRoSH5+kpdt7taWQjLbhoNCQl5/fq1v7+/tbW1yO9BBQLBDz/80LJly/79+9vb2wsXuyh1k+rm5mZoaMjhcKKiokQeoHqUlHDiBE6fXtyo0Ugnp58GD/6h8vsqKkJXNxnQAABoZGQY3b4tppjy7ubNmxkZvUpKfs7LW5uRERcXx/Pzg7Y207HKkZODnj1hY4OgICiWMe/vf5w4gcuXsWOHRJKJyOzZ+PABBw8ynYNRjRohOhqKinB2xsuXldpFSQk+PnjwAElJaNoU+/eLOaLcK7sQOjg4cMS2rOfFixevXr36/PnzuLg4Ozu7JUuWlNqAw+F4eXlFRUXp6uqKKUP1cLmf+PzQT5+upaefj4zMffr0aeX33bRpPpc7ksOZwOUOaNBgurMzbGxw9qz4wsqd/Hzs2oWZM00LCx8BPCBVR6dYT6+i8sKcDx/Qqxfs7BAYCG5FvdbevsW0aTh4UHqLepkUFbF7N375BWlpTEdhlKoqdu/GlClo1w6nKj0tuakp9u/HgQNYswbduiEmRpwR5Vu5//9CQkImTZrUp0+fHv9V81MeOXJk+PDhOjo6AH788ccjR458u82IESMaNmwovmJcbV9+YfH5VetA269fv8TEO4GB7V+/Do+Lc0tMRMuWGDwYxsbYtk30OeVKfDx+/RX16+PIESxe3GLhwm4mJm2srQceOLCe6WjlSk9Hly7o3Blbt1a8lHxJCYYPx/z5aNtWIuFEqmVLjB6N2bOZziEFvL0RGopffoGDwwIDg5ampq1OnTpf4V4uLrh/H8OHo0cP+PggNxfx8fGHDh2Kj4+XQGY5UXYhnDVrlru7e2hoqDhWH3z9+rWNjY3wtY2NTXZ2dm5ubrWP9uHDh/v/iBHzRyYNDQ0vrz4GBt00NfuWlOhaWNhWaXcTE5NJ4sL5vQAAIABJREFUkyZZWFgAMDHBqVP48AH9+8PHB9ra+O23z0ugsSpJIMCVK+jXD506AcCdOwgLg6cn/PxmJSX99fLlzS5dOjGdsWypqejeHb17IyCgUtvPnQt9fUyfLuZYYrNkCW7dwvmKf+fXfo6O2LHjrydP/s7IeJScHPnTT79VZi9FRXh7IzYWAGxsrtvbj5kwIaV9e68rV66JN67cKKPhiM/nb9u2berUqRs2bFCo1uJaMTExG8rqN71ixQp9ff1Pnz7VqfN5mJSqqiqAjx8/alV3wdNHjx5N+meVNlNT0zLvL0Vo8eJZU6eOKi4u9vdvMGECb8+ewhoecM0arFwJf3/lP/5QWr2aM2JESUBAUZ3vjSKTnMLCQkVFRcUKH16JARF950NYWhrn4EGlnTuV9PVp/PiS3btLhCMCqzS+hSkpKZx+/VQHD+b9+mtxZQJfuqR48qTKjRv5eXniHcXL+4c4Dr5+vcLUqSpOTgXq6rVkLPL3L9HvyM9/r6ZWPycHgCaPp/jhw4dK/hdTVsayZbhzZ8/t24GAXVFRL3//5e3atalGhprIz88XCARS2FxXJoFAoKKiUuFmZfwDZGRkFBQUTJgwoXpVEICOjk67du2+fV8YyMjIKPufUaZZWVlcLtfQ0LB6JwLg4uJyULLP4jU0NAAEBqJTJ+zapeHjI4Jj+vvD3x9btmDJEqWgIKV+/bB9O+rWFcGRa0LxH5I/NRFxOBzhj/pr9+9j/XqEhGDIEJw7h5YtOYAKUPGFLiXevoWbG6ZMwaxZyoByhdsnJmLaNJw8iXr11MWdTVgF64jnU1jfvnBxQUCA+po14jg8A8q7RCvk6upqbOwvEMwvKXmvqNhBVVWnEr+o/2VnZ3z37nM+3w54xuWaVCNADXE4HDU1NRkqhHw+v+Ltvu1IyuPxjIyMzp49K9purF/Mnz9/xIgRwteHDx+2s7Mrb0s9Pb0bN25851CSHD7xrTdvyNiYwsNFfNhDh8ja+vOEpcyOwZKe4RO5uRQYSHZ21KgR+fv/Z2k3GZKQQNbWtH59ZbcvKSFnZwoIEGem/5xO9MMnvpaRQSYmdPeu+M4gUTUZ6lpQUBAcHBwRcWPIEBo+vLIzsQllZGQ4OvY2NHRo2rS3pWX6jBkVDz8VrVo5fKLscYSHDh2ys7NL+P4k6tWVkJCgpaW1Y8eOq1evWltb79y5U/h+v379Dhw4IHx97NixwMBAdXX1uXPnBgYGlrd8JbOFkIjCwsjUlMSxeuO1a9Sy5eeVDp2dJ3K5Zlyu2YABEl3qjalCeO7cOVVVGwUF8969hz97Rr6+pK9P7u4UFibDI9KePSNzc9q8uQq7+PqSm5vkvmVxF0IiCgoiOztJ/+IWE5FMh1tURK6uNHVqNXfPyKA+fahjR7H8CipPrSyEZXeWOX36dEpKSuPGjVu1aiXyXqP169e/cOHCxYsXV65c6evrO3HiROH7Tk5OZmZmwtexsbH3798fNWqUsC9MYWFNH8WJiasrJk+GpydEPm9a1654+BC3b4PHS4qKihEIXgsEr8+du/nmzRsRn0n6DBkyo6AgmM9/dfFiXvv2UWpqePIEwcFwda24g6V0iouDqyt+/x1Tp1Z2l4sXcegQ9u+X1W+5TKNHo149rF3LdA6poayMkydx6xZWrKjO7vr6OH8eQ4bAyQmXLok6nDwp9/GPnZ2d+M7q7Ozs7Oxc6s3ffvu3A5Wfn5/4zi5av/2Ghw/h64s//hD9wZ2csHnzi27dtAQCRQACQZMOHfIHDcKMGWjYUPSnYxaPh7t3cfUqSkoANAUAdB4+/PD//lf6UpEtDx7AzQ2rV2P06MrukpqKiRMRFMT8c2KR27QJTk4YPBgNGjAdRTpoaeHiRXTqBEND/HNTUAUcDnx80KoVRo7EmDFYtgzV7doh3yRwcyo+jDeNCmVnk40N7dsnloPz+XxNTRsOZxKHM1FDo6m3N5mbE4dDamrUrh1t3lyd6U8rSTJNoy9fUmAgeXqSjg5ZW5O3N+nptQZmAwe43HoPHjwQdwCxunePjI3p1Kkq7MLnU/futHSp2DKVQwJNo0Jr1lCXLjLcyi0k2pXC4uPJ1JROnqz+EVJTydWVunal5GRRhSqbHDWNsqpERwenTmHOHLFM/cDlcjMynq5YYRMQ0CQ7+3FgIN69Q34+Vq0CgDlzoKICCwuMGYPHj0V/djFJScHx45g8GfXqoWNHXLkCV1fExODlSwQG4v37qEmTClxcDkdGHm7VqhXTYavvxg24uSEwEAMHVmr7hQtXmps7Ght3y819PG+emMMx5+efUVSEvXuZziFNhPNM/fgjqr3Aj6EhLl5E585o0wZSMzGl7CivQt66dWvYsGHNmze3tbUVvrNu3bpdu3aJrFKLgpTcEQodPEgNGzIwv/OjR/Sd28Rdu3YNHDj0ZLU+atbkjjA2Nnbz5q3R0dFf3snNpbAw8vUlBwcyMCBPTwoMpPImdRf3CvUSEB5ORkZU+YVGoqKi9PQ8AR7wwtq6ozijlU1id4RE9OgRGRlRSopkziYW4rhEr10jQ0OqYSPI+fOfl1sR021brbwjLLsQnjt3TlFRsUWLFh4eHmZmZsI3t27damFhIVU/AqkqhET000/k4cFYm09BAW3eTO3akZoacblkbk62tis4HGdgD4fTYvXqdVU/YDULYUTEDX39dhzOTl3dXnPnHvb1JWdn0tYmV1fy96d79yruLy6jhZDP558+fXrnzt3Hj2fVrUtXr1Zh30OHDqmorBauUG9sXMEK9eIgyUJIRL6+NHy4xM4memK6RA8fJnNzqmGH/bdvqX178vCgrCzRpPqaHBVCa2vr0aNH8/n869evfymEz58/B5AoyY66FZG2QlhcTB070ooVTOcgCg+nYcOIw+kJvAYIuKukNKhvX5o8mVatojNnKh6hGBQUpKlZX1/f5m5FI794PEpKoocP6cIF2ruXVq6kJk1mADcAApI0Nd39/Cgysmqd5mW0ELq7e2lqzlJQWKeg4HDlSm6V9n3zJkVJqRWHs1NDw2fEiOp2qK8BCRfCwkJq0oTENlxZ7MR3iW7ZQg0a1PR2uaSEfH3J0pK+apQRDRkqhNnZ2QMGjLp2reKxq2X0Gk1LS3v16tXx48e5XO7X0wcIxzakpKR8GeTAKkVJCcePw8kJ9vbo1YvJJJ07o3Nn3Lv38eXL24AlcFNLSz8tDXFxOH4ceXkoLgYRFBWhqgpNTdStCwMD1K8PCwvY2kJJKW7MmDnAGiDbyalfSkpyWhqSk5GairQ0JCUhLQ2pqUhORno6srJQty4MDWFqCkNDGBnB0tL85ct7JSXOXO7dPn0sFi1i8kchMXw+/86dZx8/RgNQV0/j8W4CVbgI/PyM+vUL6djxlKWly4AB/cUWU1qoqGDbNowdi65doanJdBppMmUKEhPRrx+uXUO1541RVIS/P9q1g4cH5s2DSCbAkjkNGjhnZo7s3t20a9cKtiyjEAqLH1HpKQGTkpIAqKmpiSZjLWVsjKAgDBuG6GjUr89wmNDQfU5OHp8++erpaT19eq1UX/wPHxATg9hYvHyJhASkpCA8HDk5+PgRxcV7ADdgJADgrK1tupmZgYkJjIxgaAgzM7RqBWNjGBvD0BCGhqXHuhUUTBs8eNLDh45WVmYbNsjUAno1QKSQn88DUgH9OnUempmNqPy+q1bh4UPcuGGqpjZNfAmljYsLunXD3Lm5Xl6xTZo0kbZl1xj0++9IS8OAAQgJQZUmYCtlwAA0a4YhQ3D3LgIDoS72efqkyI0byMoiYMH795VYzaDM+0QrKysfHx8iCg8P/9I06uPjY2BgwOPxRHfnWlPS1jT6RUAAtWlDhYVM56iu48ePA42ATCABMPP3LyoqkmgAmWsazcggV1dq0yaifv12pqb2S5f+Ufl9L18mU1N6+1Z86SpFwk2jQrdvxyootNLWnmZk1Pr2bVG34omTuC9RHo8GD6YRI6o2AVuZCgrI25uaNKEnT0QQTJqbRlNTadEiatmSlJWJyyUOpzXwcMOGimtW2YVw//79AEaNGuXn52doaHjixAlPT08A6ys/T6JESG0hFAjI05MmT2Y6Rw20adMNMAXMhg718fSkhg3p2DHJnV22CuGjR2RtTb6+VI1Pia9ekYkJ/fmnGGJVESOFcMSIn4BwgIDHPXqMkvDZa0ICl2h+PnXqRD/9JJqj7dtH+vq0bl3alCm/Dhw4qcJn/+WRwkL45580cCDp6RFAamrk4kJBQcTnU2hoqKZmo61bQys8QrnDJ3bs2GFgYPDlxlFDQ8Pf31/avn+pLYRE9PEj2dqSlI03qZqve41euUItWlC3bvTokSROLUOF8OhR0ten/furs+/Hj9S8OW3bJupM1cJIIZwwYRaHcx4gIMrDY4KEz14TkrlEP3ygVq3I3180R3v6lOrUceVwjgG3DQ2d3r9/X42DSL4Q8vn8O3fupKamfv3m8+c0cyY1bUoKCqSkRE2b0syZZcy5WqNeo0JFRUW3bt06derUn3/+KZ2/laS5EBLRs2dkaEj37jGdo7pKDZ8oKaHAQDIyojFj6L/XpOjJRCEUCGjRIqpXr5r/xAIBDRlCU6aIOlZ1MVII3717Z2XlZGDQV0HB4dSpvyV89pqQ2CX6/j1ZWdE/axPUCI/HMzBoKxyio6bmd+7cuWocRMKFMDMzU1XVksvtzOVa/vrr0qAgcnMjLS3icEhPjwYOpMjI7+0ugkIo/aS8EBLR6dNUvz6lpzOdo1rKHEeYlUW+vmRkRP7+YnwIKv2FMCeH+venzp2r/5ngf/+jDh1Iwg9fv4ORQkhEAoEgOTl5716+nZ0U/TQqJMlL9MULMjWt2kR95bG1deFwLgPPudz2U6a8qcZ/YQkXwokTvYHfAQI+AG01NKhTJ9q6lSp5qVayEJY96XZERERxWespaGlp1a9fvybr6MqbAQMQFYURI3DxYi2ZDFdXF/7+mDABs2dj9278/js8PZnOJHF//40BA9CpE44fh5JSdY5w9iz27cOdO1CueHXeWo7D4RgbG48di7NnsXw5ZGe+fclp0ABnz6JPHxgYoGPHGh3q4sWg2bOXpaZmTp++7PDheg4O2LMHTk4iCipqt2/j7FkLQBUAoKigkPrxo3jOVGZ5NDIy+s4uzs7O8cyuGPsP6b8jJCIej3r0oEWLmM5RdRXOLBMWRs2bU/fu9PixiE8tzXeE58+TkRHt3l39Izx9SkZGdOeO6DKJAlN3hF8kJcnSowTJX6JXr5KhIT18KMpjHjtGxsbk61uF1h3J3BEGBZG1NXE41KxZhpJSfS53EJfbeMqU2VU9To2aRk+dOqWvrz916tTQ0NC7d++eO3du5MiRZmZmZ8+e3bJli6mpadOmTaVhHIVMFEIiSk0lCwvaty/9+vXrqeJ+vCY6lZliTfjg0MSEvL0pLU1kp5bOQigQkL8/WVjUaLaOrCxq0ID27BFZKlFhvBAS0d691LKlbDSQMnKJHj5MFhb0119p165dSxPR/7fUVBoyhGxtK3tVi7UQ5uXRzJmkqUlKSuTmRq9fExEVFBQcP348Li6uGgesUSFs3759QEBAqTenTJkirDp3794FEPn9Z5QSISuFkIj273/I5dpra882NHSIirrJdJxKqfxco5mZNGMGGRiQv79ofotJYSH8+JEGDSJn5xotc8PjkZsbzZkjuliiIw2FkIgGDCA/P6ZDVAJTl+jcuXcVFFpra88xNBTlyMtjx8jIiLy9KS+vgi3FVAjfvqWBA0lBgbS0aOZMkX0Yqv4yTFlZWbdu3fLw8Cj1voeHR0hICABHR0cTE5OEhATRNdDWfsHBgQLBppyc1Wlp+xYt2sR0HBHT08P69YiMREQE7Oxw/jz98cfWTp2GLFiwosyHzTInPh7t2qFuXVy7BmPj6h9n7lzweNVcjlxObN6MLVtw/z7TOaTVs2fb+PzAnJyAtLRdS5duFdVhPT3x+DEKCmBnh/BwUR21Ui5fRqtWsLTEo0c4dgw5OVi7VtLPzssohEQEID4+vtT7wueCwtfKysp16tQRd7jaRFtbncPJAABkaGnVzpmOGjdGSAg2bcLkyQfmzv3rxo3Va9cWLVjgz3Sumrp0CR06wNsbgYE1+v954ABOn8ahQ7Wk25SYmJrC3x8TJ6KkhOkoUklbW53DSQcApCUkqIuw84ihIfbvx7p1GD0akycjL09kRy4Tj4eAABgbw80NOjp49gwvX2LQIPGetFxl3ie2bdvW0tIyPDxc+FeBQHD27FkdHZ3hw4cTUVpamoKCQrUnJhAhGWoaTUpKatiwnYFBVyUlp1WrEpiOUynVXoZp6NCfgNsAAWn29m7VOIL0NI2uW0dGRnT9ek2P89dfZGgomjmuxERKmkaF+venxYuZDvFdTF2iiYmJDRq0MzTsUr9+u1Gj3pmZUWBgdaY0+o6sLPL2Jmvrsi/76jWNJiUlOTn1MjJyWLlyVUoKjR5NKiqkpkajR1N2ds0jl6tGzwifP39ubW0NQFNT09raWlVVFUDr1q2FHT2ioqKmTJnCdpaphpycnGfPyMBAxF2/xKTahXD//sNaWqOBpwoKcw0Mlleja4k0FMKCAvLyInt7evOmpodKSaF69ej4cVHEEhupKoTv35OBAd2/z3SO8jF7iebmfl7k69496tiRWreuYFx5NZw/T+bm5O1Npb7L6hXCunVbAFuAv4B2HM5fRka0cqXIon5HjcYRNmrU6MmTJydOnHj8+HFKSoqFhYWDg8OAAQMUFRUBdOjQoUOHDpK8ba01tLS0tLSwfj2GDsXdu9DSYjqQeIwZMzw/v/DUqeVdujhaWf3k4YH+/bF2rQxMfl9cXPzzz4uuX7/Zpo1zbKxfkybKUVFQVa3RMUtKMGwYJk7EkCEiSikHvjSQ3rlTzZGatZvmPytXOTggMhLBwRg9Gi1aYONGkS1607cvnjyBry9atsSOHejWrToHefECp08jMhIZGQRMAQD80L37prCwXaJJKSoSqMniI3N3hF9MmkTDhjEdoiLVviMsRdjSYmND165VdhemPm4vXOhfp84SII/DWdKrl2g+snp704ABJGXT9JZBqu4IhdzcaOlSpkOUQxoaLb6Wl0f+/mRoSL6+pe/haujCBbKw+PfW8Pt3hH/9RQsXkosLmZuTigpxOKSlRU2bkpJSR+Aw8IrDaRcUFCTKfN9V/V6jLAnYsAEvXmCryPp8STVdXQQGYsMGjBuHyZMhrrkhROHq1SeFhUMBNaJhfP6jmh9w61ZERmLfvtJLNrIqY+dObNqEJ0+YziEL1NTg64v795GUhKZNsX8/vllStprc3D7/E9jZQVnZRkOjCZdrNnToaACFhTh7FpMno317GBtDQQGOjtiyBR8+oE8fHD2K4mLk5ODpU8TG7m7YcJOWlvusWV1Hjx4tmmQi9KUkHjlyxMTEZOPGjURkZ2dnUg4x1u6qk907QiJ68ULa59EQ1R3hFx8+kLc31a9Ply5VsKXkP27HxNCYMaSltU9FZQQQoa09fMeOaq0o8ZUbN8jYmF68EElAsZPCO0Ii2rmTWrWiSnymlzRpuyP8WkQE2dtTp04ifsjaocNIYDBAQDJgp6pKAKmoUP365O5OK1fS39I3a3qVnxFaWVl5eno2adIEgIeHR25uLnPVWS40aIAdOzBkCO7dg74+02kkQlsbgYGIiMAPP6BVK2zdKhXfeEwMAgJw+TJ+/BHv3nmFhqpfuHCqb1/PoUNr1JX77Vt4emL3bjRoIKqk8mjiRJw8iYAALFjAdBTZ0akT7t3D7t1wd0ffvvj9d3x30swKvHmDmzdx+zYePYoDBgMAjIHCrVvRrx/09EQUmlEcEtX9MxMOHDgQGhp68OBBpoNU38yZiI/HuXPS2HRWWFioqKgo7CElWvn5WLIEQUHYsAGDB5exARHl5eVpaGiI/NRfe/IES5ciPByzZmH6dKipiezIhYXo3BkjRmDmTJEdU9x4PB6Px5PC8cHv38PBAWFhaNGC6ShfkcwlWkN5eVi1Cps2Ydo0zJsHFZVK7VVSgsePceMG7t9HZCQKC+HoCAcHaGpe/+WXkcAM4J6KyoPCwldiji8CAoGAz+crVdThin1GyLCAAGRnY80apnNIlpoa/P1x/DgWLsTQoUhPl3SAx48xdCh69ICDA16/hq+vaKogEV24cOHgwYNeXh8bN5alKijNzMywdCnGjmWH2FeZujr8/HD7Np4+RfPmOH4cRHTt2rVjx47l/XfAfFISgoPx66/o2BG6uvDywtOncHXF+fNITkZwMPz8MHt216iok82aHRw8mCMTVbDyyr0jPHv27OrVq2NjY9XU1BITEwEEBAQoKyv//PPPkk34PbXgjhDAu3do2xZHj6JTJ6aj/Jf47gi/OgX8/bF1K5Yuhbf3v++L7+P2o0dYtgwREZg5EzNm1HRoRCmDB/9w5Uqd/HxzLvdUYuJ1AwOpHy/yFam9IwRAhD590Lkz5s1jOso/ZOKO8GtXr2LmTKSmTiko4PN49UxMgvftu3b/vvr9+7hxAwUFn2/7OnZEhw7f+1yYl5enpqbGkcImrLJU8o6w7OETe/fuBeDq6jpu3DgzMzPhm/v27TM0NJSGcfRfyHRnma9duUJmZjWazVkcRN5ZpjwPH1Lr1tSnD719+/kdcfREePCAPD2pXj1at47y80V7bCIiHo9naNhGuPy3pua8ixcviv4c4iSdnWW+SEyUrql5pLmzTHkKC3kaGk7CSxSY36hR6PTpdPAgJSRU4SASXpi3hqo/fIKI5s+f7+PjExYWNm7cuC/vOzs7p6WlvX//XoTlmiXUvTvGjcOoUeDzmY7ChJYtER2Nzp3h6Ijt20GEmJiYkJAQHo8nkuPfuoV+/dC/P5yd8fw5fHxEfCMopKCgIBCUAJmAQFk5xsTERPTnkGPCBlIvL7aBtPpUVBS0tEqAbID09J4dOWK0YQNGjhTZGHzZVUYhTE1NTUpKGj9+fKn3jY2NAaSlpUkil/xZsgQKCli2jOkcDFFUhK8vrl7Frl0wNPyhZcvhI0cGaWs3rEbv5ZcvX/7yi9/Klevz8vJu3kS/fhgxAq6un0ug+Fr+oqPB4602MnIzNnb48UdnOzs7cZ1JXk2aBENDuXugLlo7dwaYm/cwNLSfMKGVvb0903GkRRmPf5SVlQEUFBSUev/169cAtLW1xZ9KHnG5CAqCoyPatUPPnkynYUjz5rh5E0pKl4niiZTz8+cPHLh/2LBpAHR0wOFAWxtcLjQ1oagIdXUoK0NNDSoqqFPn801eZmZmp05DU1IWKCq+W7FipL7+2fnzceqU2KfpevwY/fvj4MFuffrcEe+Z5BiHg8BAODqiXz80a8Z0Gtnk5tbj3bseTKeQOmUUQj09PVtb2y1btrRt2/bLE1EiWrlypbm5eQN2VJTYGBnhwAGMGIG7d2FmxnQahigogMsF/3MbcXFBgYpwabrsbADIyYFAgNxc8Pn49AklJcjPR1ERCgpQWAgVFSgpRefl9SUaVFICoiMxMTxVVTF29hF68QJ9+mDDBvTpI+5TyTtLy88NpLdvs3OQskSm7N8R/v7+AwYMeP/+va2tbUFBwcaNG48fPx4ZGRkUFCQrnYVklIsLpk2DpyfCw+X3//nYsYP27m0BmGpoJF+58qjyAxuKihAT06hnz4CsrA9AYt26kEAVfPsWPXti8WIMHSruU7EAYPJknD6NtWvh68t0FFatUV4vmuDgYOEsM0Lm5uaSnCm1kmpNr9GvCQTUvz/Nnct0Dgn2Gv3Wu3fvrl69Wr19Dxw41rhxp/btPWJjY0Wb6luJiWRtTVu2iPs8kiDlvUa/9vo11a1LMTFMZpDFXqMiUSt7jVYws0xiYmJqaqqWllaDBg2k8F6wdowj/FZ2NhwcsGYNBg5kMoYExhGWh2RhkFZ6Orp0wdixmDuX6SiiIM3jCL+1dSv27UNUFBQUmAkgE5eoONTKcYQV/I4zNzc3NzcXXarP4uPj16xZk5mZ2atXrwkTJpT6mX78+PHYsWNRUVFFRUVt2rTx9vZWFUdvdymmq4sjR+DhgVatYGXFdBpWWXJy4OYGT89aUgVlzo8/4swZrF2LOXOYjsKSfQxMsfbx48dOnTrp6emNGjVq1apV69evL7VBVFTU2bNnO3ToMGDAgIMHD44cOVLyIRnXpg3mzcOwYSgqYjoK6xv5+ejXD87O8PNjOoq8EvYgXbw40NLSuWtXT2Gfdharehho9Tp48KCNjc2yZcsA1KlTZ/LkyTNmzOBy/y3JvXr16t27t/B1s2bNWrRokZ+frybCGZFlhI8PbtzA7NnYtInpKKyvFBdj0CA0aIB165iOIt+Sk28Rhb59e+Xt2weDB/94//5FphOxZBUDd4R3797t9M+smp06dXrz5k1KSsrXG3zdUvru3TsdHR15axr9Ys8eXL2KAweYzsH6R0kJBg+GpiZ27JDGBUPkSlzcs+LinoAq0CE5OZPpOCwZJpY7wpKSkuTk5G/fNzY2VlZWTklJadq0qfAdNTU1NTW1lJQUU1PTb7fPzc2dMWPGkiVLvvNg9ubNm926dRO+NjExCQwMFMV3IEX27uX26ZOwc+cqLS2On9/0Ro0aSezUzHaWyc/Pl/x5v4/Pxw8/1CHCjh2F30w4IfN4/2A6SGW1bdtGV3dsero5h/PA3LzFp0+fJHl26bxEJSA/P18gEMhQZxmVSqw+JZbfcc+fP/fw8Pj2/ZMnT9rb26uqqhYWFgrfEQgERUVF6uplTNKfn5/v7u7erVu3qVOnfudcjRo1mvPP43JVVdXa14mrVasigWBMePhmQPDXX5MSEm5XPJO6iCj+QzKn+xoRcTgcqfrXJMLkyfjwAcHBqFNHioKJimz1GgXQrFmza9f27tp1TEXFfPfueW/eKEtyuhkpvEQlg8PhyFyv0Qo3E8vvuObNm796Ve5qVRYWFm/evBHG+MMqAAAgAElEQVS+TkxMJKJvbwcLCgo8PDxsbGw2b978/Z943bp1XV1da55ZaiUkJCgr2wLtAfB4jV+/ft2wYUOmQ8mjX35BTAwuXxbjbKWsqmrevPkffzQH0KwZBg3CnTtgp4BkVQMDzwg9PT2Dg4MzMzMB7N27t3fv3pqamgDCwsIePHgAoLi42NPTU09Pb8eOHV93opFPVlZWiopPgUggoqjoeX12ongmzJ+Pa9cQEgL5uwGQDWPGoEcPjBkDgYDpKCwZxECZ6dChQ//+/Vu1auXi4rJ9+3Z/f3/h+2vXrj19+jSA06dPh4SE3Llzp3HjxjY2NjY2Nl/uIOWQiorK1auHBg0K6tPnkJLS4Vu35HXiNeYsX47gYFy5Al1dpqOwyvfHH8jNxdKlTOdgySAGHv8ACAwMfPXqVUpKSuvWrb88kzh8+LDw6dfAgQOzsrK+3l7Ol7xo0qTJyZPbAVy9ipEjER0tv1NyS96mTdizBxER0NdnOgrru5SUcPQonJzQvDkGD2Y6DUumMFMIAVhbW1tbW3/9jo6OjvCFsrKycCkoVindu+OnnzBkCP78E5XoCcWqqX37sGoVwsPBLrIrE4yMcOYM3NzQtClsbZlOw5Id8v4ETub8+issLDBzJtM55MCpU5g3D5cuset3y5LWrbFqFQYNQk4O01FYsoMthDKGw8Hu3YiIwK5dTEepjQQCwerVm3v2HO3tvXnqVEFICL5agoUlG7y80LUrxo7FdxcUYLH+xRZC2aOhgVOnsGAB7t5lOkqts3bt1sWLn4WFLdy5M2748EB7e6YDsaplwwZkZ2PZMqZzsGQEWwhlUqNG2L4dnp5IT2c6Su1y6dLNT5+mA02IZsTE3GA6DqualJRw7Bi2b8f580xHYckCthDKKg8PjBmDYcMgO1NiyQBV1fYczhYgQV19c+/eHZiOw6o+IyMcO4bx4xEXx3QUltRjC6EMW7wYqqqYN4/pHLVCURHGjcO7d1PnzrXq3Nl34ULrWbOmMB2KVSPt2mH1agwahNxcpqOwpBtjwydYNcfl4sABtGkDR0cMG8Z0GlmWkYHBg2FkhKgorpqaD+DDdCKWaIwdi5s3MXYsTp1iVwthlYu9I5Rturo4dQozZiAmhukoMuvxYzg5oX17HD0K+Vv1svbbtAmZmVixgukcLCnGFkKZ16IF1q2Dhwcy2RXZqi4kBD16YNky+Puzdwy1k7DjzLZtuHCB6SgsacUWwtpgxAgMGICRI1GJ9UZY/1q/HlOmICQEI0cyHYUlTsbGOHoUEyYgPp7pKCypxBbCWiIgAMXF7IzDlVVUhLFjsW8foqLg6Mh0Gpb4tW+PRYvg4cF2nGGVgS2EtYSiIk6cwP79OHWK6ShSLzkZLi4oKkJUFCwsmE7DkpQpU+DsjHHj2BlnWKWxhbD20NfHyZOYMoUdOPU9jx6hQwf07o3Dh6GqynQalmRt2YL0dAQEMJ2DJWXYQlir2Ntj9WoMHMi2/5Tt5En07Ik1a+Dnx3aNkUfCjjObNiE0lOkoLGnCFsLaZswYuLiwMw6XRoSVKzFzJi5cwKBBTKdhMcfEBEePYvx4vHzJdBSW1GALYS20aRMyMrBqFdM5pIawa8yZM7hzBw4OTKdhMa1DByxciAEDPm3YsHPnzt35+flMJ2IxjC2EtZCSEo4fx8aNuHiR6ShSIDkZnTujuBjXrsHYmOk0LOkwZQo/MbH3zJm5P/2U2bZtH2LbT+QbWwhrJ+HAqfHj8fIlZcrxSPtHj9C+Pdzc2K4xrP94+fKlsnI9gWBWcfGc1FSDN2/eMJ2IxSR2rtFaq0MHTJz4wtZ2hI6OvrZ2bkTEaWP5uCFKSkravfuQoaGepuaon39W2bYNAwcynYklZYyMjIheAJ8AQWbmy8uXDSZNYvtPyS/2jrA2e/x4RXHxhrS0S/HxvosXr2M6jiRkZ2e3aeP+v/8ZTZuWNGHCmAsX2CrIKoO2tvaGDfPr1+9uZdVr5crF+/apt2mD6GimY7EYwt4R1mZ5eQWADgAinY8fC5iOIwm3b9/OyelLNKakBNra7Vu25LEXOatMw4cPHD7886ek2bMRFISBA+HqijVrYGDAbDSWpLF3hLXZkiXTDQ29dHRmKijMtLefzHQcsSsuRliYdX7+baAAeKWuzldUZKsgq2IcDry88PQpTE3RogXWr2en7ZUvbCGszZydO8TGXjx3btD163/6+9veucN0IHG6fh329njxovHixaPq1+9iZzf55MltTIdiyRIdHfj7IzwcoaFwdMSNG0wHYkkK+3m5lqtbt26nTp0AbN+OQYNw61YtnF0zKQm//oqoKKxfD3d3AOMWLhzHdCiWrGrcGBcvIjgYY8bAyQlr1tTC/zKsUtg7QnkxcCCmT0f//sjLYzqK6PB4WL8ednbQ1cXjx8IqyGKJQL9+iI2FrS3s7eHnh6IipgP9Izs729m5v4mJk7OzR3Z2NtNxagm2EMoRX184OWHYsFry/CM8HPb2CA7GjRtYvx7q6kwHYtUuamrw88Pt27h3D3Z20jI9xZw5K6Ojh6Wk3I2OHjl79gqm49QSbNOofNm0Cb17Y/58rFzJdJQaSE6Gry+uX8eyZfDyYjoNq1Zr0ADnz+PKFUyfjs2bsXEj6tev6THv3bu3bFmgkZHekiWzDQ0Nv90gPx9v3iAxEe/f4+1bJCYiMfHzi0+fUvj80QD4/Jb794ecPQsLC9Srh3r1YGHx72tTU5TZUezjx4/r129PTc2aPn1so0aNavqd1BZsIZQvwtnX2reHtTUmy2A3Uh4Pmzdj6VKMGoW4OGhoMB2IJR9cXfHoEbZuRdu2mDIFs2YV//rr7Pj4d0uW/NquXbsqHSo5Oblv36lpaRu53Hfh4SN37LiSnIxXr5CUhC8vsrNhagoTE5iawtoazZqhRw9YW8PEBI8fjx49+scPH0bo6BwJCvqtXbv/7Bgejlev8OoV3r6FpiasrT/vJTyOiQlmzBgZG+tWXNz+xIkRjx9fNGBHigBgC6Ec0tNDaCg6dkSDBujenek0VREZiWnToK+PiAjY2jKdhiVnlJXh44MhQzBvHvT1O/B4HYGBV68OvX//XKtWrfLykJv7+U9ODj58QG4uPn78983s7M8vkpIepqf3BtoKBG1fvFj7228l9eopWVigaVP06IF69WBmBn39cmOYmrpGRJhGR0e3bRtoa2sLQFcXzZqV3ozHQ1IS3r7F27d49w6JiYiKwps3/JiYdIFgKoCPH+9HR0e7s8/VAbCFUD5ZW+PIEQwbhvBwyETrSGYm5s9HSAiWL2fbQllMMjPD/v04eDATWAdAIMhq1+4Rj9dKRQVaWp//6OhAW/vfvxoZoWFD6Oh8/mtBQYthw5ZkZo7icN42aqR0/bpSVTPY2traVvRJUFHxcxvpfyk0aMB9+fIhUK+g4HpQ0GgnJxgZVfX8tRBbCOVU585Yvhxubrh9W0rn0cjPz9fQ0BAIcOAAfH0xdCji4qCpyXQsFgtQVKTi4meADYcT5u8/ZNq0sh/IlcP8+HH/xYt/MzLSX736gPhClun8+d3e3gsyM7Nmz54RF2fTvDnmzoWPD5SVJRxEunBkev2RAwcOhIaGHjx4kOkgsmrOHERHIywMKiplfLWwsFBRUVHyk7MsXRrg57cBUK1b16R+/XBVVc6mTWjeXMIp5A6Px+PxeHXq1GE6iAw4c+bsyJGziouLu3Rpf+XKMabjVF98PObPx4MHWL4cnp6V2iUvL09NTY0jIzOUCwQCPp+vpFTBbTdbCOWaQIDBg6Gjgz17yvgqU4VQSakej/cI0AWmDB9uf/iwt4QDyCe2EFYJEeXl5WnUiv5a167h559hYIA//oCdXQUb18pCyI4jlGtcLg4dQmwsAgKYjgIASE/Hpk3g8TiANgAOx8TE5BnToVisWq5bN/z1F4YNQ8+e8PJCWhrTgSSOLYTyTlUVZ85g0yYcPcpYhsJCBAdj6FA0aICICDRpYsfldudwZioqbp8zZw5jsVgsuaGoCG9vxMXB1BR2dli5EsXFTGeSIMYK4fPnz589K/fDfl5eXkxMzNOnT4vl6l+DIaamOHsW06ZBwrNy8/m4cQOTJ8PEBCtXwtUViYk4dgxxccFHjkxduFAzJeWxiYmJRDOxWHJMV/fztOORkWjRAiEhTAeSFAZ6jebl5fXt2/f9+/ccDsfExCQkJKRUO/v169cHDhxoZWVVVFSUnZ198ODBbt26ST6nXLG3l+is3LGxCArCvn3Q04OXF+LiYGz8nw2GDBni5uZWOx7AsFiypXHjz5Pp+Phg3TqsW1fGOMVahoE7wsDAQIFAEBcXFxcXx+Vyt20rvVZOmzZt0tLSHjx48PTp05kzZ/r4+Eg+pBySwKzcb95g5Uo0boy+fQEgPByxsfD1LV0FWSwW41xd8fAh3N3RpQt8fPDhA9OBxImBQnj06NHx48crKioqKChMmDDh6DfPptTV1ZX/GdVibW3Nrx1TRMsCEc7KvXHjDlfXkUuWrC4pKcnKwvbt6NgRDg549Qq7diEhAf7+sjGWn8WSW0pK8PGB8BFW48ZYvx69enloajarU8fq0qVLVT3au3fvrl69yuPxRB+0xhhoGn3z5o21tbXwtbW19du3b7/dpqioaNGiRTk5OXfv3t28efN3jpaRkXHlyhXha1VVVWdnZ5EHlisimZV7z56DCxfezM31i4zcFxQUkJa2oFcv+Pqid29U1I2ZxWJJF319rF+PH3+Eh8fm+PgE4Hxxcayb20SBILHyB5k+3XfLllOAtYrKxMTEv/T09MQXuBrEUggvX7584EDpGRMUFBT27NkDID8//8tYpTp16nz69KnMg+jq6ioqKubk5Dx+/Lhr167lnevvv/9evny58LWxsbFdhaNgWBXZs4fTvbuaqWmxo+N9XV1di4qeGfL5SE7mvH3LffuW+/o159077sWLN3NzfwQaFRf/XFLi9eLFJ1VVAlBYiMLCSmUgovz8fJke5CpzhOMIS0pKmA4iG+TtEjU3h4rKVmAc0BxoTrSKywWHAy4XAJSUwOEQl/v5k66KChQUwOWSmhoAqKrSo0chRA8B9YKC37t0WTN69G/m5lS/vsDKSlBhTwBf3/8dPHjZxETv1KnACn8dlfLp0ydVVVUdHZ3vbyaWQmhpadm7d+9Sb3K5n5thjYyMsrKyhK+zsrKMy3pApKKi4uvrC2Do0KGOjo4TJkzQLGdyrQ4dOrAD6kVLUxOXLsHGxh4QAIVGRjopKY+EX8rO/jy3/Zc/SUl4/RoqKp/nube2RseO0NPrsHPntrw8fVXVPcOHdzQ0rHKfFyLicrlsZxlJYgfUV4m8XaI5OVBXHwXsBXoBsUDGy5f4+BHCtYGTksDncwoKPv81IwM8HoqLOTk5AJCbi0ePBP8cifvypeHy5XUKC1FSAoEAHA4UFKCiAlVVaGhAWxt6ejAwgJERzM2RkHAoMPAeUXBu7g0Xl2Hp6U8qn3n69F+3bDm0atWBWbM6f39LsRTCxo0bN27cuLyvOjo63rx5083NDUBUVJSjo+N3DqWlpcXj8QQCwXe2YYncgQNLAGPgOsBPTbVr0SIvO1s9PR26urCygo0N6teHkxOGDEH9+rC0LD1R4adP/UNCNr5966ajozZ1ajBD3wSLxRKN4GD89BPc3Obp6ESHhfVTUkJIyE4rqyocYerUPoGBrYBGysqxb97cq1v33y+9e4f4eCQk4N07JCcjNRUZGYiNxa1b+PQJubnviH4ArAHrjIwNwttQBQVwOFBUBIcDFRVwOFBVBQANDXA40NCAggLU1HD9+mmiZ3x+WRNI/hcDzwinT5/ep0+fZs2acTic9evXBwd//kXZqFGjXbt2derUaffu3fn5+Q0aNPjw4cPatWsHDx6sra0t+ZzyLC0tDRA24isAynz+R319dQ0NfPqEp09x+za0taGp+fmPlhZ0df/9q6YmwsM3vn49nsebnJJyeebM30+e3M7w98NisaolNRVz5uDWLezbh65dAZyp3hRrW7as9vWdHhsb27Nnz1KzNgrXEy7v8dehQxajR/9B1BoINzAQREcjJwc5OcjLQ1YWiouRkQEipKYC+Pz6wwfw+fj4EURcoE5BQcXxGCiEHTp0OHLkyK5du4jo0KFDHTt2FL7v4eEhXKy5WbNm+/btCw0N1dbWHj9+/IQJEyQfUs6NG7d682YrwBPIUlLKfvq0dPP1hw/4+PHfP8K11oSvMzLw4kUqj9cRAFGzsLBNy5ahd2/Y24PLTmTEYsmO48cxfTrGjcOTJ6h5k7mlpaWlpWVV9xo5cuTffyds3uxlZWUUHBxapUWjunRpGhnZXkMjEGj1/S3ZSbdZpcXHo2tXrFjBu3bN29jY+EtfpMq7c+du375Tc3LGaGqenjp1RknJwCtX8PYtunSBqyvc3WFqWsERatOMxrKCfUZYJbX7Ek1IwOTJSE/Hzp1wcPjPl2Rr0u2wsLDGjRvX+2ZhxlLY9QhZ/5GYiJ498b//YfRoxSFDtlRv6Yk2bZzu3DkWGRnZuvXG5v+sn5SSgsuXcf48fH1haop+/eDqChcXdkAFiyVFBALs3ImFCzF7Nn75BQoKTAeqme7du1dmJDpbCFn/Sk9Hjx6YOhWTJtX0UFZWVlb/fZhubAwvL3h5gc/Hw4cIDsavv+Lvv+Hign794Ob2n6nd3r17Fx8fz06tx2JJ0pMn+OEHqKggMhLl93eshdiHNqzPPnxAr14YNQq//CLeEykowMEBfn64dw8vX8LLC/fvo21b2NjAxwdXrmDs2BlWVt179FispdUwPz9fvGlYLBZQUvJ54vuJExEeLl9VEGwhZAnl5cHdHZ07Y+FCiZ7XwACenggMxLt3OHQIenr47Tfs339GIHgiEIR//Oi+YsUKiQZiseTPzZto1Qo3buD+fXh7Q0Ye/4kSWwhZKCiAuzsaN8YffzCWQUEBbdti0SLcuiXsXCq8MtVSUthHiCyWuOTn49df4ekJX18EB8PcnOlADGELobwrKcHQoTA3x44d0vJJcMiQXlxuKy7XXVk5+Ny5uT4+KGcaPhaLVWURERFnz57Nz8+/eBG2tnj1Co8fw8uL6ViMYjvLyDU+H15eUFDAnj1SNMjv6NEdsbGxcXFxAwYMyMxUnDMHLVti0ya4uTGdjMWScWPGTA8J+VRcbKGg4K+vf3XHDrUePZjOJAXYQii/iDBlCjIyEByMao2SECNbW1tLS0tFRUUjI+zfjz//xJQp2LIFW7ZIYt1gFqtW4vP5ly7dzs6+C0BFhbduXWSPHr2YDiUVpOYugCVxc+YgJganT4tgzghx69IFDx7AwQGOjli/XgTLJbJY8oYIZ88qZGfzgByANDTizczqVrybfGALoZxasABXriAkBLIyM0adOvDzQ1QUzp+HkxPu3GE6EIslO27dgosL/vc/zJq13NS0i4FBq1GjmjiUmjNGjklZixhLItatw4kTiIiAri7TUaqoQQOEheH4cXh4oH9/rFoFLS2mM7FYUuzpU/j5IToaCxZg4kQoKLitXMk+bC+NvSOUO1u2YNMmXL+OKk1fK1U8PREXhzp1YGuL/fuZTsNiSaV37zB5Mrp2hYMDnj+Ht7fMz5cmPmwhlC8HDsDfH2FhFU97LeV0dbF+PY4cQUAA3N3x+jXTgVgsqZGZiV9/hb09dHXx99/w9ZWBfgDMYguhHDlzBnPn4uJFVGlFTWnWsSMePECPHmjTBn5+KC5mOhCLxaj8fKxciSZNkJ2NmBj4+4Ndy7Uy2EIoL65cweTJOH8etrZMRxEpJSX4+CA6GnfuwNERN2+ioKAgODj4xo0bTEdjsSSnpATbt6NhQ9y/j1u3EBgI49KriLLKxXaWqeVOnz4XHn7Pyqrz8uWuJ06gdWumA4mHlRUuXEBwMIYPL8jJ6S4QdFVWTnR3P75v33qmo7FY4kWEEyewYAEsLXHuXOnlA1mVwRbC2mzjxp2//XYtJ2cMh7N24cKCTp36MZ1IvPr1A59/Y9iwTsXFywBcvNiGx+NVb0lFFktqZWRkxMbGNmvWrG7dulFR8PVFURG2bkX37kwnk1ns74ja7PDhkJycjYA5kVFMTCBQywshgHr19DU1X2VmEvAxP7+Yw2GvcFatcv/+/b59vYuLu3A4s1q3DnzxwnH+fPzwgxRNkSiL2B9ebWZo2JTDOQvwVVXPOTo2ZTqOJLRu3XrEiMYGBq2MjFwaNFju6oq0NKYzsViis3jx1tTUbdnZa7KydqSmbnnxAt7ebBWsKfbnV2udOoXbtxd07/7UyqrD6NF5c+b8xHQiCdm48fe0tEcpKQ/u3evj4gInJ0RHM52JxaoZPh9XrsDbG2FhakAWAA4nq1kzVSV2mTJRYBuOaqdNm+Dvj0uX1Fu23Mx0FsYoKMDPDw4O6N8fS5di0iSmA7FYVSQQ4OZNHD+OY8egpwdPT1y+PGfcuCF5eavV1XMDAk4wHbCWYAthbUOExYtx9Chu3kS9ekynkQL9+iEiAoMHIzISgYFQVWU6EItVkW/rX2QkGjQQftHi5cvo7OxsXZmbIFGKyXbTaHR09B129uWv8Pn48UdcuICICBFUwfnz5+/evVsUuaosLCxs7Nixojpao0a4dQtFRejYEQkJojpqbbNnz5558+YxnUJmXL16dcyYMaI9pkCAGzfg4wNzc0yeDF1dREQgNhZ+fl+q4GcMVsFx48ZdunSJqbNXVWhoaGV+k8j2HSGPx+PxeEynkBZFRRj9//buPaiJc20A+LtJuAQCuXETgwGVixhFpirUGEACqAgCar9KpaV4l3Z0dNpRy/QMI61QqJVWrX6ohXC0o/RAsTelSkwEolaORYOIfkC8UFFBCEIgQJL9/thpTg6CIqIL5vn9lfeSNw9J2Ce7++6+CaitDZWWIju7ERhQq9VqtdoRGOj59fX1aTSaERyQwUDHj6OcHDRnDpJIUETECI79miDx4x6LXuQrqtVqz5w5w2Kx5s6di0z2//71L8RiobfeQjIZ8vIa0XBHjkaj6evrIzuKoert7R3KxzS2EyEwUqvR4sXIzQ2dPIng/Plg1q5FU6ag+Hi0ejX6xz9grh0gQXd3t79/2L17wRYWf82aVeDt/Y0x/0mlyNub7PjMEiTC18H9+2jhQiQSoexs2Lg/g0iELl5Eb72F/v1v9M9/IhaL7ICAeejtRXfvotu30alTFbdvC7XanQih0tIAoVAnk9E8PcmOz7xhOI6THcPwpaWl7d27d/r06WQHQrLa2k10ehOfXzCyw6pUKjqd7kLGLQvVavX9+/d9fHxe0vgGA62ubi2G6T09//clvcSY8+DBg66uLo/X5o7sL1N3d/fly3d6erDJk20nTBj/lJ46nW19/ftarUtPDwchjErtolIvPH78fxTKxzSaBsP+JyBgjM15uXHjhrOzM2uM/IRsa2vz8/M7fPjw07uN7USIEFIoFF1dXWRHAQAAYDQKCAiwe9akiTGfCAEAAIAXASeUAAAAmDVIhAAAAMwaJEIAAABmDRIhAAAAs0ZNTU0lO4bhaGhouHr16u2/tbS0uLq6kh3Ua6K5ufnSpUvu7u7GmrKyMktLy2fOvBqKe/fuXb9+ffz4/8w4l8lkLBbL2tratNuJEyeuXLkiEAhMK/V6vVQqdXV1Jdba1Wq1MpnM2tp6RAIzZzdv3qyurjb+N6nValKumRkr7t69e/fuXScnJ7IDeUWqq6sfP37M4XCIYllZGY7jxssnzp49y2az+/3/ko7YVnC5XGNg165da25uHvRTw8emLVu2jBs3LvRvq1atIjui14dWq/X19c3LyyOKxcXF48ePb2trG5HBDx8+PHv2bNMae3v78vLyft3S09NTUlL6VXZ0dCCEbt26heN4Z2dneHh4dHR0d3f3iARmzlatWuXm5mb8b9q4cSPZEY1q6enpUVFRZEfx6qSlpS1evJh4rFaraTTau+++SxQbGxsxDHv48CF50Q0qKSkpISGBeNzY2Mjlcp/czhiN4TvLzJ8/Pzc3l+woXkNWVlaHDx9etGhRaGgonU5ft26dRCIZVdfPtrW1LVq0aOLEiXl5ecTeIXhBS5Ysyc7OJjsKMBqFhITs2rVLr9dTqdSysrKIiIiysjKiSSaTCQQCR0dHciMc0O7du6dNm1ZUVBQXF7d69eo1a9YIhcLBOsNGBAwgMDAwKSlp/fr1NjY2sbGx8+fPJzui/3jw4EFMTIxQKNyzZw8FbigHwEs2e/bsvr4+pVI5Y8YMuVweHR3d1NTU0NAwceJEuVweEhJCdoADYzKZOTk5K1euvHnz5p07d4qLi5/SGbYjYGBpaWnXrl37448/srKyyI7lv0RHR8+bN2/fvn2QBQF4BSwtLd98802ZTIYQkslkwcHBQUFBcrkcITSaEyFCaMGCBfPmzUtJScnNzbWysnpKT9iUgIHV1ta2t7d3dna2t7eTHct/mTVrVklJyYMHD8gOBABzERwcLJfLOzo6/vrrLx8fHyIRNjU11dXVBQUFkR3doDQazcWLF21tbRsbG5/eExIhGIBOp1u9evXOnTuTkpLWrVs3giPb2dkRc14Ier2+q6vruaZ97t27NyQkJDg4uKmpaQQDAwAMJiQk5Ny5c3K5XCgUYhgWFBR09uxZuVwuEAgcHBzIjm5Q27Ztmzp1alFR0YYNG5qbm5/SExIhGMCOHTuYTOb69es///xzlUqVn58/UiP7+vo2NDQY9+cuXbpkYWExud/y20+FYdi+fftCQ0PnzZsHuRCAV2D27Nk9PT179+4NDg5GCDk4ODAYjLy8vNF8XLSioqKgoCAnJycsLGzRokVbtmx5SmeYLAP6q6qq2rNnT2VlJYZhxAzSqKgosVhsevHfsE2dOjU+Pj48PPy9997r7e3dv39/amqqjY3Ncw1C5MLk5OTQ0FCpVDpu3LgXDwyAoVMqlRs2bCAeM5nMjIwMcuN52YjThCUlJS+WkNIAAAivSURBVJmZmURNUFDQgQMHRvZw0QjSaDTvv//+119/7ezsjBDavXu3QCAoLi6OjY0dsP9YvaCexWJNmzaNz+eTHchrqKam5p133vH39yeKbm5uU6ZMoVKpI5VvYmNj+Xz+7du3qVTqRx99tHz58if7lJeX9/b2hoaGmlZiGMbj8QICAiwtLTEMi4yMZLPZOp3O9Np/MAxsNtvPz4/H45EdyNhgZ2c3adKkcX/j8XjmsCSqp6dnQEDAggULMAxDCHl4ePj5+UVHR4+2S+kJKpXK29s7Pj6eKFpbWwcHB2s0Gm9v7wH7wzJMYDTKyMjo7Oz87LPPyA4EAPD6g3OEAAAAzBrsEYLR6OHDh3q9Hk7+AQBeAUiEAAAAzBocGgUAAGDWIBECAAAwa5AIAQAAmDVIhAAAAMwaJEIAwDMUFRWdPn2a7CgAeFlg1igA4Blmzpzp4eHxww8/kB0IAC8F7BEC8Oq0tLRotdrhPVetVo/U4lPNzc2tra2Dtba2tj569GiIQ7W1tT18+PDpHdra2p4vPgBeLUiEAAzTyZMnORzOzZs3ieKhQ4c4HI7x5r0dHR2Ojo4HDx5ECD169CguLs7Ozs7R0dHW1tbHx6ewsJDo1tDQwOFwcnNzTUcuKiricDh//vknUTxx4oSvry+bzXZxcZkwYcKxY8cGjOf333/ncDjnzp0zrczMzHR2djamvQMHDkyYMMHJyYnL5QoEAmK1VaODBw9OnDiRy+U6ODhwudwvv/wSITRjxoyqqqqffvqJw+FwOJyEhASic0lJydSpUzkcjrOzs7u7+/fff28cJzs7m8PhXLp0yc/Pj8PhLF269PneWQBeMRwAMCytra1UKvXbb78lisuWLbO0tJw1axZR/OWXXxBCSqUSx/Fbt26tWbPmt99+q6mpKS8vX7JkCY1Gu3z5MtEzMDBQKBSajhwZGenl5UU8/vHHHykUSmJi4vnz5y9fvpycnIxh2KlTp56Mp6+vz8XFJSkpybTS29s7JiaGeJyZmUmhULZv315ZWXnhwoWYmBg6nV5TU0O0ZmVlIYRWrFghk8mUSuXRo0ezsrJwHC8rK/P09BSJRKdPnz59+nRVVRWO4xcvXrSwsJg7d65MJlMoFLGxsRiGFRYWEkPt3LkTIeTh4bFr1y6FQiGVSl/wrQbgpYJECMDwzZw5c+nSpTiO6/V6BweH5ORkKpXa2tqK4/jmzZudnZ0NBsOTz+rr6+PxeFu3biWK+/fvRwjV1tYSxfv379NotIyMDBzHDQbDpEmTIiIiTJ8eFBQUFhY2YDybN29mMBgdHR1EsaKiAiFUVFSE47harba1tf3www+NnbVaLZ/P37BhA47j7e3tDAYjKipqwGHfeOONZcuWmdbExMSwWKz29naiqNPpPD09p0+fThSJRLh///7B3jcARhU4NArA8InFYqlUqtfrq6qqHj16tH37dltbW+J4Y2lpqVgsJtasQQh1dHQcPHhw69at69at++CDDwwGQ319PdEUHx9Pp9OPHj1KFI8cOWIwGFasWIEQamhoqK+v9/LyOmOCx+MplcoB41m5cmVnZ2dxcTFRlEgkXC43MjISIVRRUaHRaHg8nnGcsrIyPp9fXV2NELpw4UJnZ+eqVauG+IdXVVUtXLjQ3t6eKFKp1GXLlimVStPTgTExMUN/JwEgESzMC8DwicXiL774oqqqSiqVCgQCHo8nEolKS0tFIlF1dfWmTZuIbleuXBGLxXQ6PSwsjMvl0mg0CwuL9vZ2opXJZMbExOTn56emplIolPz8/Pnz5xOrAxKzY/Ly8oxp0shgMFAo/X/ICgQCf39/iUSSkJCg1WoLCgoSExOtrKyMQ6Wnp/d7loeHB0KopaUFITTEJQkNBkNjY2O/W6K7urriON7a2spms4kaYk1UAEY/SIQADJ9IJLK2ti4tLZVKpWKxGCEkFotzcnJEIpHBYDAuLJydnW1vb69UKm1tbYmaX3/91XScxMTEY8eOyeVyJpN59erVTz75hKhnMpkIoV27dq1du3aIISUmJm7ZsuXOnTsKhUKtVicmJpoOVVxcHBIS8uSzWCwW+jtZPhOFQqHT6UTuNGpubja+irHbEGMGgFzwTQVg+KytrQMDA0+ePFleXm5MhLW1tRKJZNKkSe7u7kQ3lUrl4+NjzIJ1dXU3btwwHSciIsLNzU0ikUgkEiaTuXjxYqLex8fHycnpuS7gS0hIsLCwOHLkiEQiIXYQifo5c+ZYWFgUFBQM+KyAgAArK6vBWhkMRnd3t2lNYGDgmTNnenp6jDU///zz5MmTHRwchh4qAKMF2ScpARjb0tLSEEI0Go2YOWIwGJycnBBCa9asMfbZtGkTnU4/deqUVqutrKz09/dnMBjh4eGm42zbto3BYDg4OKxfv960PicnByG0cuXK2trarq6u+vp6iURCTKUZTFxcHI/Ho1KpX331lWn9xx9/TKFQPv30U5VK1dXVVVtb+8033+Tm5hKtW7dupVAoKSkpKpWqs7NToVBIJBKiKTk5mc1mFxYWVlZW1tXV4TheUlKCYdjbb79969ate/fubdy4ESGUk5ND9Ccmywzr7QSABPBlBeCFKBQKhJDp9Q/Lly9HCB0/ftxY09LSIhQKiZ+eVlZWaWlpYWFh/RKhcR/x/Pnz/V7i0KFDLi4uxh+vjo6O/TJcP8RkGRqN1tTUZFqv1+t37NhhnOGCEOLz+QUFBUSrTqdLSUmh0+lEE41G27ZtG9HU2NgYGRlJHD6NjY0lKr/77jvj6UAbG5v09HTjC0EiBGML3GINgFcBx3GVStXa2url5WWaiobIYDDcuHGjo6PD2dmZ2NsbdiR9fX3Xr1/v6elxdXUdP358v1ZiT5FCofD5fGOeG0xvb29NTY1Op/P19bWxsRl2SACQCxIhAAAAswaTZQAAAJg1SIQAAADMGiRCAAAAZg0SIQAAALMGiRAAAIBZg0QIAADArP0/WoxuY7grHiwAAAAASUVORK5CYII=", "text/html": [ "\n", "\n", "\n", " \n", " \n", " \n", "\n", "\n", "\n", " \n", " \n", " \n", "\n", "\n", "\n", " \n", " \n", " \n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n" ], "image/svg+xml": [ "\n", "\n", "\n", " \n", " \n", " \n", "\n", "\n", "\n", " \n", " \n", " \n", "\n", "\n", "\n", " \n", " \n", " \n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n" ] }, "metadata": {}, "execution_count": 7 } ], "cell_type": "code", "source": [ "plot_bandstructure(scfres; kline_density=10)" ], "metadata": {}, "execution_count": 7 }, { "cell_type": "markdown", "source": [ "or get the cartesian forces (in Hartree / Bohr)" ], "metadata": {} }, { "outputs": [ { "output_type": "execute_result", "data": { "text/plain": "2-element Vector{StaticArraysCore.SVector{3, Float64}}:\n [-4.79982738197534e-16, -8.864207065081506e-17, -3.337615924791053e-18]\n [-1.690636556339103e-16, -2.575110771448328e-16, -6.890003712560586e-16]" }, "metadata": {}, "execution_count": 8 } ], "cell_type": "code", "source": [ "compute_forces_cart(scfres)" ], "metadata": {}, "execution_count": 8 }, { "cell_type": "markdown", "source": [ "As expected, they are numerically zero in this highly symmetric configuration." ], "metadata": {} } ], "nbformat_minor": 3, "metadata": { "language_info": { "file_extension": ".jl", "mimetype": "application/julia", "name": "julia", "version": "1.8.1" }, "kernelspec": { "name": "julia-1.8", "display_name": "Julia 1.8.1", "language": "julia" } }, "nbformat": 4 }