{ "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", "[density-functional theory](https://docs.dftk.org/stable/guide/density_functional_theory/)\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": [ { "output_type": "execute_result", "data": { "text/plain": "3Ɨ3 Matrix{Unitful.Quantity{Float64, š‹, Unitful.FreeUnits{(ā„«,), š‹, nothing}}}:\n 0.0 ā„« 2.7155 ā„« 2.7155 ā„«\n 2.7155 ā„« 0.0 ā„« 2.7155 ā„«\n 2.7155 ā„« 2.7155 ā„« 0.0 ā„«" }, "metadata": {}, "execution_count": 1 } ], "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.900404478290 -0.70 4.1\n", " 2 -7.905002484823 -2.34 -1.52 1.0\n", " 3 -7.905210708346 -3.68 -2.52 3.0\n", " 4 -7.905211515396 -6.09 -3.34 2.8\n", " 5 -7.905211530433 -7.82 -4.61 1.4\n", " 6 -7.905211531398 -9.02 -5.25 3.5\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.1020985 \n AtomicLocal -2.1987887\n AtomicNonlocal 1.7296105 \n Ewald -8.3979253\n PspCorrection -0.2946254\n Hartree 0.5530405 \n Xc -2.3986218\n\n total -7.905211531398" }, "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.0184081\n 0.261073 0.116914 0.00482502 0.0611642 -0.0239772 -0.0184081\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.417656\n 0.354533 0.389828 0.384601 0.449215 0.627552 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 k-points). There are 7 eigenvalues per k-point because\n", "there are 4 occupied states in the system (4 valence electrons per\n", "silicon atom, two atoms per unit cell, and paired spins), and the\n", "eigensolver gives itself some breathing room by computing some extra\n", "states (see `n_ep_extra` argument to `self_consistent_field`).\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+gvaeTAAAgAElEQVR4nOzdd2BTVdsA8OfcNN17pS3QzRJKBy3QAkIZshRQFEVAHO+LEwRFQWWLih8v0y2vKPiqDBHEhQIts6VAyygdbCmQNuluupPc8/1xtZbSQmlzc9fz+6tJT5OT05zz3DMvoZQCQgghpFSM0BlACCGEhISBECGEkKJhIEQIIaRoGAgRQggpGgZChBBCioaBECGEkKJhIEQIIaRoGAgRQggpGgZChBBCioaBECGEkKLJLRBeu3ZNr9e3MrHZbOY1M7KEhdYGWGhtgIXWBlhobSO3QLhixYotW7a0MnF1dTWvmZElLLQ2wEJrAyy0NsBCaxu5BUKEEELormAgRAghpGgYCBFCCCkaBkKEEEKKhoEQIYSQomEgRAghpGjKDYT/+27Lo/966f/WfGAymYTOC0ICyMrKeuxfL06f9bpWqxU6LwgJSaGB8IedP770wbbfus1cknT9zaXvCp0dhKytvLx8+MNTt3hO+K85PvGBR4TODkJCUmgg3LXnQPnAF6FTr+rhr2//fb/Q2UHI2k6cya4ISoAuA2nk/TrWqaCoROgcISQYhQbCxPhYl+NfQ1m+zeH1xf5xsTtN/z3HVuEQKVKAMyX0hSPmR86GmC8eBW0OXEqtM5RE7XZ+47j5ioEKnTuEBKDQQDhtyuMLx/fu/evzM8Jq8jcvWR6n+uM6DfzO+Oxh85kSbAuQDNWZYdsVdvhvplG7za62cPoJzb6vPxh0Ytn9f3555rfvDo9VU4B+u0zDfzNtu8KaWKGzi5AVEUpl1e7PnDmzc+fOM2bMaE1ig8Hg4uLS8DC/GjZdYD/OYX3sYXo3Zmo442DDW0Ylq0mhodYQttAuVdD159ivzrMRnmR6N+bBIMamhQvgOjPsymM/z2WzS2FqZ/LiPUwnJ2LdzP4Dv2ltgIXWNgrtETbL3xHmRjJXHrVZHqf6KY8N3WKch4NFSLJYCntv0In7zH1/NNWYIGWszZ5RNo+EtBgFAcBOBY+EMHtG2ewdraoxQcR20wN/mPbekNfFMkK3wC5PUwyBYR3IsA42Fyvof8+x/XaZet3pOhohUSmogY3n2U9yWG97mN6N2ThIfbdjG93dydp41bJY1XeX2Dlp5hozPN2F+Vc3xsuOnxwjJChs2lsU7kqWx6nyHlNP78Z8nssGbzHNO26+VoUXx0i8DhfQifvMPb43XjbQH+9TnRhvM71b20f4XdQwvRtz6iGbbxNVlw2081bjxH3mvTewCiC5wR7hHXCDRY+EMDll9NMcNuoH09AAZno3ZmgHItjkCUI3K6+HLZfZdVksS+HZbsyGe9XOaku+fm9v8tkA1ft9VBvPs88dMbvbwvRuzORwxgnbDyQL2CNsLW6w6M/H1MM6kDlp5m7bTO+fZovr/vqt2WwuKcGdWIh39fX15eXlDQ/Ti+izh80hW4x7b9A1/VTZD9u83JOxbBRs4G4LL/dkzj9iszxOtffGX6usMxutsi4pKcE7pCMrKCwstOwyT7yiuzvcYNH0bkx6Ef08l+281TgsgEmoSP2/N2aYnbw6OjOHftvp6OgodDaRPH2z5fvZC94mju4xnQMfW/bl2mwwGOFfXZnzj6i97a2Uh78n0VX51apNF9j7/zD72MO0wOr1L43XG20YQ+EPGz/t17evlXKDFOb69euDxkyosvNwqClO2rU1JCTEIi+L2yfatdS4sBa+PM8umDS0/qmN4O5vm7R29b2uLzw3vT2vKXK4PrsNLFVoHbrHaGckg60j+W7WgJHjFj4+XPAheiMLP15l33z/g4sGQhOfh5JrkT89e+rQnva/Mn7T2kD2hfbkC7M32Q6hPUfC+YMPF27b9uWnFnlZHnuEe/bs+emnnzw9PadPnx4QEHBrgkuXLn3xxRc1NTWPPPJIQkIC92R1dfXGjRuzs7N9fX2feOKJoKAgALh48eL27dsb/vCRRx4JDQ3lL+et52MPr/di1tsbL9o7A4DJ3q2yulroTCHZMlMKNnYA4OjiNqdr7bAOws9Tqxl4OITJDq5ffMkXAMDepbauXuhMIdmqqqmlbm4AAA5uVdU1lnpZvuYIt2zZMmXKlG7duul0uvj4eIPB0CSBVqvt27cvy7JBQUGjR4/ev38/9/yYMWN27tzZr1+/srKymJiYGzduAEBOTs4HH3xQ+jej0chTtttmyeuzvD8d47B1lsPBT5+cMkno7CDZeuSxycza+z23vhCqS73vvvuEzs4//j1tSseUdW7bZzPrHpg35xWhs4Nka/7s59WbX3be/qpmy/Qlr7Vq5K9VKD+ioqI2bdrE/TxgwIBPP/20SYKFCxc+/PDD3M8rVqwYNWoUpZRbb3Lt2rWGF/n6668ppbt27erfv39r3nfGjBnr1q1rZSYrKipamfKObty4se/A4Y4bDScKWUu9pjhZsNCUw1KF9tAe07xfLhw9etRoNFrkBS2ourr60KFDk3Zcf/O4ySIviN+0NpB9oX170RzzTeHBQ4dKSkos+LK89AgNBsOpU6eGDRvGPRw2bNjBgwebpDl06FDjBIcOHQIAV1fX4ODgY8eOAUB+fv6NGzciIiK4NHq9ftmyZR9//HFeXh4feW6ngICAIff2fyvO8Y3juGoO8eJ4IT1WSBfcF963b18bG9Etc3NwcBgwYMCK4QGf5bDXcbst4oGRhYXp7IrBXgMHDPDw8LDgK/NSnfLz8wHA29ube+jr63vgwIEmaQoKChonqKysrKiocHV13blz58iRI1955ZWioqKVK1dGRkYCgIuLS0JCAsuyhw8fnjdv3s8//3zvvfc2+9ZGo/Hbb789deoU9zA2NnbatGkt5bOmpkalUrXvs95kciCsybT55XJtop9sDy22eKEpgUUK7bWjqrciKNTXV4t4Ds6DwNRQZskJdm1ce68I8ZvWBvIutM/OM6HOpJ97/V2txLC3t2eYO3T5eAmEtra2AGA2m9VqNQDU19fb2TU9mkmtVjfcGr6+vp4QYmtrW1NT89hjj7344otTpkzJysp65plnYmJi4uLiBg8ePHjwYC7xokWLFixYcGtk5TAMExgYGBcXxz0MDw+/9a0bNJuxdlrSmy46TUYEqYVfxsAPPgpN9tpfaLuvQ0EtfaorEf85f/NjoPt29uWequ7u7aoE+E1rAxkXWqURVmSzPw1n7OzuLmyRViys5iUQ+vv7Mwxz/fr18PBwALhx48atq0Y7dOjALYThEnh5ednb2+/du7ekpGT+/PkAEBwcPGbMmG+//bYhqnH69u371VdftfTWKpVqwIABzz33XGvyqVKpLH71NDEM1mSZfrhKJoaKvsVqEz4KTfbaWWgshQUZpuV9GDu1BL5UXg7wSgRZeopuHdqu7wl+09pAxoW29gw7JIDp7cPLp+OlXtnZ2Y0cOXLLli0AUFtb++OPP44bNw4Aqqqqdu/ezXUEx44d+/3333PnUGzZsoVL4OPjU15ertPpuNc5f/68RqPh/pB7hlK6c+fOXr168ZFtiyAAy+NU89NZo2wHR5G1fXeJVTMwLkgCUZAzqydzVE+P6nGmEFlGUS18kGVeEsNXFeBryn3JkiWjRo06ffr0hQsXQkJCxowZAwB5eXmjRo0qLS11d3d/4oknNmzYcO+99/r4+Bw7doxbLBMZGfn444/36dNn2LBh2dnZ5eXl//73vwHgmWeeuXr1amBgYG5ubmVl5a+//spTti1ikD8JcYYvzrHPdZdMy4VEq56FRRns+oEqCQ2226tgQTQz77h5/xjRLepBUrTslPnxMCbMla9KwOPJMsXFxYcOHfL09BwwYAA3V1lXV5eTkxMREcF13o1G44EDB2pqagYNGuTq6trwh2fPnr1w4YKPj0+/fv241XFVVVUZGRkFBQUBAQFxcXHcHGSzrHyyTEtOl9CRv5nOT1S78HPqo4Bkf3QFH9pTaOuy2D032J/uk1hEMVPotd20Ol51X1t3/eM3rQ1kWWh/GmjsTlPWw2qNA19vgUes8fWleTzZ3MODvBUlt06hLGsa39pcaJVG6LLN+OsImygvCXUI//LDn+zSDDbjQRumTXnHb1obyLLQpu43h7uSRbyNiwLefYI/78QyqzPNeoudAYSU6D+Z5uEdGClGQQB4KJhxtIGtl3G2HLVdZgndc4N9JYLfUIWBkC8hLuSxMOb9M7i/HrVRYS18lM3yeiHMt+V9VG+eYOsxFKK2mnfcPD9axfcck4TrmPgtilFtPM9eNshq8BlZzdIM89RwJtRFkt1Bzr1+pJs7rM/FSIja4lABzS6Df3flPU5hIOSRjz28eA+zNANbAXTX/jTQ7y6x8yIlvyfs/TjVspNmg7jOyUfSMO+4+Z1Yxo7/SoCBkF+v9VL9cYM9XYKdQnR33jrBzuqp8uVtmZzVRHiSoR2Y1WfxchDdnR1/spVGeMwqJ5NgIOSXsxrm9lItOIGtALoLZ0pokpZ9uadMquey3sy6s2YdLhxDrWamMP8E+399VG1bcny3ZFLTxOz5e5isUro/HzuFqLVeP2ZeGMP7AgGrCXYhk8OZ907jwjHUWl+eZ73tYURHK02QYyDknS0Di3sz846bMRKi1jiQT8+XwzP8LxCwpgXRqm8uspcqsBKgO6s1w9IMdnkf602Qy6qyidbkMKbeDD9dxQFSdAcUYN5x83txjK28qqa3Pbx0j2oxLhxDrbAui+3rS+J9rbdeWl61TawYAstiVXOPsyZsB9Btbb/C1pjgkRAZVsxXI5h9WvZkMXYK0e2U1cPKTPPS3latAjKsb+I0uhMJcISvL2IkRC0yU1iYzv6nr5UWCFiZsxrejFLNP4Ezheh23jtlfjCIaefNLO8WBkLreS9OtTCdrTEJnQ8kVv89x2ocYFhbT6kWv2e7MefKIUmLnULUPG01/eIcOz/a2oEJA6H19PEhcT7k4xzsFKJm1JjgnZPs8jjJ76C/DTUDS3szb+DCMdSCRens9G5MRydrXwtiILSq9+KY90+bS+qEzgcSnzVZbIKG9LXiAgFBTApjTBR2/omXg6ipc+X0x6vsnF4CXAtiILSqrm5kXBDzHzyJG92stA7WnDW/bd0FAoIgAO/Gqt7AhWPoFm8cZ1+PVHnaCfDW8q94YrOkN/N5Lnu9CgeH0D/eOWV+OITp7Cbz7iBnREfS0Qm+uoCREP3jWCE9Xkhf7C5MSMJAaG0BjuSZrsyyk9gKoL/cqKKbLrALouU8O9jE8j6qxRlsNS4cQ39747h5aW/GwUaYd8dAKIA3olQ7rrI5ZdgpRAAAC9LZZ7szftI/X7v1Yr1JvC/5MBsvBxEAwK/XqLYapoYLFo8wEArA3RZejVAtTMdWAEFuGf3tGjsnQkHdQc47scx/zuDCMQQshfknzO/HMTbChSMMhMKY2YNJ09OjeuwUKt3c4+zrkSo3W6HzYXVd3MhDwcz7eBK34n17ibVVwQNBQgYjDITCsFfBwhhm3nFsBRQtTU/PlNAXBFogILhFMaovzrHXcOGYgtWzsDiDXR6nEnadmEJroBg81YUprIHfr2MroFzzjpuXxFjjBtzi5O8I07sxS/EkbgX7JJu9x50M9hd4vTQGQsGoCCyLZV4/ZmYxFCrST3lscR1MFm6BgBjMjVTtymOzceGYIlUa4f0z5rdjha8CwudAyR4MZpzVsOUyXhErDkthwQn2vTiVwENCQnOzhdd7qeafwCqgRCvOmO/rwER6Cl8HMBAKbHmc6q0TbB3OFSrMpgusiy2M6SR8EyC4F+9h0otoig47hcpSWAsfZbMLY0QRg0SRCSUb6Ee6u8P6c3hFrCD1LLwt9/O1W89eBYtx4ZjyLMkwT+vChLqI4loQA6Hwlseplp00G4xC5wNZy4dZbC9P0l8jiiZADJ7ozJTUwa/XsFOoFFcMdPMldq4Q52s3CwOh8CI8yfAOzKpM7BQqgsEI/8k0LxPBAgHxUBF4N5aZdxwXjinFWyfY2REqX9GcpoS1URTe7s18kGXW1QidD8S/90+bR3Vkenhgd/AmY4MYVzV8cwkvB+XvdAlN1rIze4go+ogoK0oW7EKmdGbePYXTJDKXXw2f5rCLxLFAQGyWx6kWpuPCMfl7Pc28KEblohY6H41ghRSL+VGqby+xlypwbEjOlp40P92VCXTG7mAzBviRCA/yaQ52CuXsQD69WAFPdxVX6BFXbpTM2x5m9FAtxlM25OtCOf3+Cjs3UiwLBETo//owy0+bK3DhmExRgHnHze/FMbYiizwiy46yvdKT2adlTxZjp1Ce3jrBvhqh8hLiBtxS0c2djOzE/OcMDo/K0/dXWBMLj4SKLu6ILkNK5qyGt6JUb53AVkCGThTRFD0V1QIBcVrWm/k0hy3AhWOyY6awKJ19T+jztZuF1VJcpndjLpTDz+crSktLhc4Lsgyz2VxQUPDGMdOiaMZRoBtwS0gHJzK1M7P0RF1BQYHQeUEWo9frPzlT09EJhnUQYRwErJfiomYgMn3tg29s9nC0HTckYf3aFULnCLVLbm7u8AmTKx39qkoLv9i/EyBA6BxJwP3sqWFTntnqG+BDK478scvT01PoHKG2M5lMw8c/mlVgKNbrli1eCDBB6Bw1A3uE4lJdXX1g13emeYcLX97/49GcS5cuCZ0j1C6vLnrv+vgPyp7daRr95rsr1wqdHWl4Y9ES9tnvip/98ULUM6s//lzo7KB22b17dwYJLnzuF3bugY9XLRc6O83DQCguJpOJqO2BMAAA9s61tbVC5wi1S21dPdg5AwC1daqqwf9mq9TV1YGtEwCwts41WGgSV1dXZ7Z1AgCwsTWbRboqHodGxcXV1fX+gX12/XdCCXWMCoB77rlH6Byhdlk27+URU56sD+nneS11/o+bhc6ONLz31pxpsyaU+Md63Eh7ec+PQmcHtcvo0aM7vr/uwtdaz7ILr7zwL6Gz0zxCKV+L9YuLi3Nzc8PCwvz8/JpNUF9fn5GR4erq2qS5v3z5sk6nCw0N1Wg0DU+azeaMjAw7O7uIiAhCWpxunTlzZufOnWfMmNGaHBoMBhcXl9Z9Gqs6d+7c6J+rt02LivEW3cSyaAtNtO7fURRVc3be2FhnZ2eh8yIZRUVFc37M9ukas2IAFtpdEGf1/Pp8/Zd7Tm4a16Fjx45C56V5fPUIt2/f/uyzz0ZFRZ06deqdd9559tlnmyS4evXqkCFD/Pz8dDpdr169tm7damNjYzabJ02adOzYsYiIiLS0tBkzZixYsAAAioqKhgwZYmtrW11d7efn98svvzg4iOa4Vh507dp1VLE5KZ+KMBCiu8JSSKtyWzMCo+Dd8fb2fnRw7yVZuOlSDg7qVA8lxnXsKN6ZOF5yZjQaZ8yYsXHjxr179/7xxx9z5swpLy9vkmbZsmVDhw49cuTImTNnsrOzd+3aBQDJycnJyclnz5796aefkpKSFi1axO0iWLlyZefOnU+cOHH69Onq6upNmzbxkW1RSfQnyVqRjqej1jtVTH3ticYeD0m4a3282exSWlYvdD5QuyVpaaK/qK/peQmER44cYVl29OjRABATExMWFvbrr782SbN169ann34aABwdHR999NFt27YBACHEwcGB6+15enqqVCqWZbnETz75JACo1eopU6ZwieUtMYA5oqNGDIUSl5RPhwSIugkQLVsG+vqSwwV4DSFteZW00kTvEfftVngZGr127VpwcHDDTF5ISEheXl7jBOXl5RUVFcHBwdzD4ODgffv2AcCQIUPGjx8/duzY6Ojo/fv3r1u3zsvLi1J648aNxomvXbvW0luzLHvu3Lm9e/dyD/38/Hr27Gnhj2cVnnYQ4kJOFNF4X1F/gdDtJWvZZ0R2vrCEJPozyfns/YF4OquE7dPSIQGMyFsxXgJhdXW1ra1tw0M7O7vq6uomCQCgIY29vX1VVRUAGAyG3NxcNzc3R0dHGxubM2fOUEpNJlN9ff2tiZtlMpmSkpKys7O5h4MGDWqIoLeqqqq6zbobwQ30sfn9T1OEo0nojNxE5IUmKiYWDhfYfhxbh4XWBlVVVf08VLNP2FT2wPPWWkuE37Q9eeoEH7ayUrCTIx0dHRnmDhejvARCPz+/4uLihofFxcX+/v6NE/j6+qpUqpKSEu7MiKKiIi7Bf//73/r6em7k85VXXunYsePjjz8+cOBALy+vkpIS7m8bEjfL1tZ2xowZrVw1SikV8xKG+wLp2izzYmd7oTNyE5EXmqik6mmoqznIy8lgYLHQ7hal9F5fx6sHjLU2zt7iqgTiJcLqeaTItLSPrbO4bz3Gy6BNTEzMxYsXCwsLAaCuru7YsWN9+vRpnEClUsXExBw6dIh7ePjwYS5BZWWlh4cH96SdnZ2Tk5PBYACAPn36HD58uEli2RvkT9L0tBaP4JasZC1OELaLDQP9NeRgAU6VS9WFckophLuKvRbw0iPs1KnThAkTpkyZ8vLLL2/atCkmJiYmJgYAvvzyy6+++urAgQMA8Morr7zyyiuurq5//vnnH3/8sXLlSgAYN27c8uXLV65cGR0d/cMPP5jN5v79+wPA7NmzH3300Q4dOlRUVHz99depqal8ZFtsXNRwjwdJ09NB4l5whVqSnM++3APnt9olMYBJ1tKHgoXOB2oTqSwW42sa/4svvhg4cOCGDRtCQ0N37tzJPdm9e/cHH3yQ+/mxxx5bu3bt1q1bc3Nz9+/fz220jIyMPHTo0JUrVz755BMXF5fU1FQ3NzcAGDZs2P/+97+ff/45LS1t9+7dyjlvJdGfJOfj5bAk1bOQpqcD/CTQCohZoj9J0uLCUalK1tJEKQRCHk+WEYQ8TpZp8McNuuyk+eD9IjoJT/yFJhIH8unc4+ajY20AC61NuEJjKWi+MWZOUPvJ+QgNixHVN40C+H9jTBtnEyTuCULAQ7dFboCGnCymVeJaN4paJTmfFfkmYklgCAzQMHi4hBRllVJnNRF/FAQMhCLnaAPRXiRFJ6teu0IkaWliANYvC0gMIMn5WAWkR0KLxbCiit2QADxrTXqqTXCymPbXSKMVELkhASQZpwklKDlf7CerNcBAKHaJ/kwSXg5LzREdjfYiTiKa25WwHh7EYKRXK7EWSAlL4WA+K5UV7xgIxS5eQ/DoYclJ1rJSGRQSPwIw2J/Zj5eDknKqmGocSICjNGoBBkKxs2Wgjw8ePSwxSfk00R8rl8Uk4uio1CTlS2PjBAfrqgQkBjC4m1BCDEbILqV98bR0y8HdhJKTrJXSqmkMhBKAiwWk5UA+7eND7PFIGcvp4kYowMUKrAXSYGIhRUcHSWdQRDIZVbI4b3LZQItqhc4Hap3kfBY3Tlhcoj9eDkrG8SIa7EIkdFQ6VlcJwKOHpUVC26ckBHcTSojkqgAGQmngjh4WOhfozkrq4FIFjfWWUisgCUMDyD4ti3VAEpLzWWktFpNSXpUMFwtIxf58doAfUWPFsrRAZ+JsQ3LKsBaInRSPm8f6Kg3RXkRfSwvwTt2il4wnq/EmMQAvByUgVUfv8SDutkLn425gjZUGPHpYKpK0kjlWSnJwvYwkSPG4eQyEkoGLBcRPXwP5NTTKS2KtgFQMCWAO5LM4TyhyUjxuXmLZVTLcTSh++7TsID9GhXGQH/6O4ONATpdgLRCvahOckuBx8xgIJQOPHha/ZEkdKyVFeDkocod1NEqCx81jIJQMPHpY/CS3fUpyEv0JHjcoZhI9bh4DoZTg0cNidq2KVhhpDw/ptQISkhjAHCqgJgyFYpUszePmpZdjJcPdhGKWpKWJ/gyGQV552UGwM0kvwlogRtI9bh4DoZTg0cNilqzFCUJrSAwgeKtqcTqQT/v6SvK4eQyEEoOdQtHan487CK0h0Z/ghlpxktzJag0kmWklw92E4nSxghpZ6OKGgZB3g/yZo3paZxY6H+gWSZIdFMFAKDFDA0gSHj0sPslaOlSaTYDkuNlCN3dyrBArgbiU1MFlyR43j4FQYvDoYXHCHYTWNAQPHRUfSR83L81cKxu2AmJDAZK10jtfUboS/RncTSg2kj5uXqr5VjLcTSg2OWXU0YYEu2AgtJKBfiS9iFabhM4HaiRJS4dI9loQA6H0JPoz+/HoYTGR7hoBiXK0gUhPkqLDOiAW3HHzkZI9bh4DofT4O4IvHj0sJsl46yWrGxKAZ62JiNSPm8dAKEl49LB4sBQOFrCDMRBaV2IAg1VAPKS+WAwDoSTh0cPicbqEetuTDk4SbgWkKMGXnC2lFUah84EAQPrHzWMglCTu6GEjhkIRkHoTIFF2KojzIYcLsFMoPBkcN4+BUJK4o4cz8OhhEUjOx40Twkj0Z/CsNTGQwXHzGAilCo8eFgMzhSM6Okia5ytKHVYBkZDBcfNYgaUKjx4WgxOFNNCJ+NgLnQ9F6utDLlXQ4jqh86F4MjhuHgOhVOHRw2KQJPHFcpJmw0C8LzmIq8YEJY/j5jEQShV39HAaHj0sKDxZTViJAQzejEVYSbI4bh4DoYThbkJh1bOQpqf34gShcLAKCE7qOwg5WIclDI8eFlaannZ1J+62QudDwaK9iLaa6mqEzodSUYD9shgUwUAoYXj0sLCScAeh0FQEBvox+/FyUCCyOW6er0BYV1f38ssvh4WFxcbG/vTTT82m2bRpU2RkZOfOnRcsWMCyLABcv359+M1+/fVXADh69GjjJ48dO8ZTtqXF0QaivPDoYcEk57OJOC4qtER/HB0VjGyuBW14et1ly5ZlZGTs378/MzPzscceO336dEhISOMEx44dmzVr1q5du/z9/cePH+/v7//CCy94enrOnTuXS3D9+vVnnnlm/fr1AFBYWJifn79mzRruV6GhoTxlW3K4s9aGdVAJnRHFqTVDehFN0MihFZC0xADyUQ72CIWRrKUTQuRQBXi5nqWUrl+/fsmSJZ06dRo9evTo0aO//PLLJmk+//zzJ554YsCAAWFhYW+++eZnn30GAI6OjsP+dv369SFDhgQHB3Pp3d3dG37l7QgLjVwAACAASURBVO3NR7alKDGAwZv0CuJwAY30JC5qofOheBGepKKeXq/CWmBt3HHzg/wwELagqKhIp9P17t2bexgTE5Odnd0kTVZWVkxMTEOCnJwcbnSUQyn96quvnn766YZncnNz+/XrN3r06I0bN1KKX/q/JPiSLDx6WAjJ+aw8BoWkjgDc64ebKAQgp+PmeRkaLSoqIoS4uLhwD93d3fV6fZM0xcXFrq6u3M9ubm5Go7G8vNzDw4N7Jjk5ubi4ePz48dzD8PDwTz75JCwsLDc3d86cOaWlpbNmzWr2rWtra+fOnbtw4ULu4fjx4z/44IOW8llVVUWI5P+LMZ7qPVeqRwRYaXRIHoXWfnuvqRdFmisrW1XsWGht0PpCS/BS/XGVPOiHy8as+k3bfUU10IdUVtZa5+3azNHRkWHu0OXjJRB6eHhQSisrK7lQV1FR4eXldWsag8HA/WwwGFQqVUNcBIANGzZMmTLFwcGBe9i9e/fu3bsDQExMTF1d3Zo1a1oKhPb29osWLZo+fTr30NbW1snJqaV8UkqdnZ3b+CFFY1hH9mipzYQuVpomlEehtVOlEbIrjImBdg6tq0BYaG3Q+kIbFUJX55qdnfGkO6t+01JKTNM6M/Iodl6GRn18fFxdXXNzc7mHOTk5ty5v4bp3DQlCQkJUqr+a8vLy8h07djQeF23M3d29puZ2+4YcHR09/nabKCgbQ/DoYas7WED7+JBWRkHEt65uxEzhsgFrgfWYWDhcQAfLZdU0Lx9DpVJNnjx5xYoVJpPpwoULO3bseOKJJwBAp9O99NJLXBibNm3apk2b8vPzuR4el4DzzTffhIeHR0dHNzxz6NChoqIiALh69ep77703cuRIPrItUX3w6GGrS9bixglxGexPcNWYNaUX0SBn4i2H3iAAf/sIly1bVllZ6e3t3a9fv0WLFkVGRgKAwWDYuXOn0WgEgBEjRjz55JPdunXTaDS+vr5z5sxp+Nvk5OSZM2c2frU9e/aEhoY6OTlFR0fHxMS8++67PGVbimwYSNDg0cNWhWdtiw3uJrSypHyZ7CDkEF5XYNbV1dnZ2d0mAcuyJpPJ1rZVp1TV1tba29/hCmTmzJmdO3eeMWNGa17QYDA0rOiRtBVn2GtVdF28NaYJZVNobVZWD4HfGYumqm1bfRmJhdYGd1VoVwy0/0+mG4+r5dM2t4nVvmn3/Waa0YN5IFAm4yL8fozbR0EAYBimlVEQAO4YBRUrMQDHhaxnfz7bX0NaHwWRFYS4EDsVOVeGtcAauOPmB/rJpw7I55MoWbQXya+mBXj0sFUka2liAFYc0RmC04TWclR2x81jfZYD7ujhAzhNaBWyOV9RZhIDCG6rt45k2VUBDIQygYsFrENfA9eraLSXrFoBeRgSQJK1LIuVgH/yO25eVh9GyXA3oXUk57P3+jMqjIPiE+BIvOxJZinWAn5VmyC9iA6QxRGjDTAQykRPT1JRT/MqsRXgV7KWyuA2pHKFN6y3ghQdjfIiTvI6TQIDoUxwRw/vx04hz5LltX1KZhL9cZqQd8n5crglfRMYCOUDFwvwTVtNi2tpTw+5tQKykRjAHMxnzVgJ+JQkx1XTcvs8SjYEdxPybJ+WDglgGIyDYuVjD52cSUYR1gK+GIyQVUr7+citDmAglA/u6OFLFdgK8CVZiyeriV0i7ibk08ECGifH4+YxEMoKzpHwKjkfV8qIXWIAScYNtbyR63HzMvxISpaIq+Z4c9lA68y0mzsGQlEb7M+k6mg9hkJ+yHWxGAZCWRkaQPZpcUsxL5K0dIjs1gjIj7stdHYjx/RYCSyvrB4uVdA42U0QAgZCmQlyJo42JBePHuYB7iCUiiG4fJofyVo2QUPUcgwacvxMyoaLBXhyoECeg0Lyk+jPJGtxbNTykvOpLCcIAQOh/OBuQj7kllEVgRAXDIQSMNCPHC+iNSah8yE7SfJdNY2BUG6GBjD78ehhS0vS0mEybQLkx1kNvTxJKk4TWpS+Bm7I97h5DIRy4+8I3nj0sKUl58v2WliWEv1xE4WFJeezg+R73DwGQhnCI2YsiwIcyGcH40oZ6UgMYLAKWJa8F4thIJQh3E1oWWdKqIcd6egk21ZAfvprSGYJNRiFzoeMJMl0ByEHA6EMDQlgDhXg0cMWI7/7ccuevQp6e5MjOqwDlqGtpqV1tId8j5vHQChDXnbQyZmk49HDFoInq0lRYgBuorCYvTdoor+cj5vHQChPeIdSSzFTOFTADpLp9ikZG4Ibai1H9ovFsHrLE66as5SMItrRiWgchM4Hukt9fMn5clpaJ3Q+ZGG/rCcIAQOhXA32Z1J0tM4sdD6kL0nWi+VkzJaBfr7kYAFeDrbXpQpaZ6Zd3eRcCzAQypObLXR1I8cLcWiovZLzWXkPCslYYgCDpyy1X3I+HSr34+Zl/vGULDGAJGEr0D5GFlJ1dJAfVhNJwplyi1DC/aixhssWHj3cfml62tmNeNgJnQ/UJr29ybUqqqsROh8St18Bq6YxEMrWQD+SjkcPt49cb0OqECoCAzQMThO2R04ZVTPyP24eA6FsOashwpOk4NHD7ZCsZeV63xmFwFOW2ilJS4cq4FoQK7mcDQkgODraZrVmOFFE+/vJvxWQMbw9ZzvJfgchBwOhnCX646q5tkvR0QhP4qoWOh+oHSK9SEkdvVGFtaAtKMDBfFb2E4SAgVDeEvDo4XZIVkYTIG8E4F68HGyrMyXU0450UMBx8xgI5cxeBbHe5HABtgJtkaSliXLfPqUEif44TdhGSYo5bh7rucwlBjB41lobVJkgs4Qm+CqiFZC3IQFkHwbCNlHCDkIOBkKZs8va/emsSdOen6XX64XOi2Rs3ra93+iHXX58o66yTOi8oPZyLM/T//fF/g88dvjIEaHzIhnnz58fN/mZ3Uue8Cg4JXRerOEOgbCqqur8+fPl5eXWyQ2yrJycnP/7v+WG+5f8jxkwdvIzQmdHGtLS0l5c/vnZxHcLPHs+9q8Xhc4Oaq/7JkyuiXwwJe6NCdNnFxQUCJ0dCWBZdtiDk3Z1esI4bPbUp/9lMBiEzhHvWgyEX331VZcuXZydnbt27eru7h4YGLhmzRpKcYRBSjIyMgw9x4FvOBs97uqNfKGzIw2pacdLIx8FnxC23+SzOeeFzg5ql7q6urI6Ct2HQsA99Z0TMzMzhc6RBOh0unoXfwhPgE6RpqC4c+fOCZ0j3tk0++yKFStef/31vn37rlixQqPRFBcX//LLL7Nnz87Ly1u1apWVs4jaLC4uzvW9J4vC7yXarPCgTkJnRxoG9o/32PBaSVC0+uLh6F49hM4Oahc7OztvR5vC0z9RV43d+X2Rka8KnSMJ8PPzs6/SQdYesHNUXz3RrdtKoXPEP3qLuro6V1fXZ599tsnzb7/9tkql0uv1t/6JeMyYMWPdunWtTFxRUcFrZsQgKXl/3/HTAh55o6SkxCIvqIRCW/6/X90HTZ0zf4nBYLDICyqh0CzOUoWm1WqfeukV9aCn/0hJt8gLipmlCu3KlSv2w18cM2V6VlaWRV5Q5JrpERYWFlZUVLz4YtPZkRdeeGHBggVXr1718fFpTXzdtm1bRkZG9+7dJ0+ebGPTzBulp6f/+OOPzs7O06ZN02g0AJCbm3vw4MHGaSZOnOju7g4A2dnZ27ZtU6vVU6ZMCQwMbH2kV7jEwYP29B8U8K3RyQ13hrcW7TniqeX3reinEjojyAL8/f03fLCy4HdTFW6GaTWTZ5DPU2t+ntT8kKH8NPPNcHZ2JoQUFxc3eb6oqEilUgUFBbXmdefMmbNs2TKNRvPZZ59Nmzbt1gRJSUlDhw61s7O7ePFibGxsaWkp9xbpf9u2bdurr77KRdCMjIz4+HiWZQsLC3v37n39+vW7/qAK5qKGMBdyqhjnd1srVU/jNYpYNa4c8b5Mqg6rQGul6GiCoqpAs/3EYcOG9ezZ8+LFiw3P3LhxY9CgQTNnzmxNN7OoqMjBwYH785KSEgcHhwsXLtz6FqtWreJ+vu+++1auXNkkwfTp05966inu50mTJr3xxhvcz5MnT274+VY4NNqs5w+b1mSaLfJSsi80llKfr+uvVbIWfE3ZFxofLFtoe2+wA34yWvAFxclShfac5VoMSWh+rODLL7+sr6/v0qVLTEzM6NGj+/TpExISkp6ertVqJ/4tKSmppeB69OjRjh07hoWFAYCHh0dcXFxycnLjBCzLHjhwYMSIEdzDESNG7N+/v3GCmpqarVu3PvPMXyv+k5OTR44c2ZC4yauhO4rXkFS8DUXrXCynDjakowKOlVKUfr7kVDGtx7MlWidVYT3C5oeAO3bsmJ6e/vXXX+/fv//69essy0ZERADAlStXGtJUVFS09KL5+fm+vr4NDzUaTX7+TWv3i4qKjEZjQxqNRqPVahsn2Lp1q6+vb0JCAgCYTCa9Xt8wMXlr4saMRuO33357+vRp7mFsbOyTTz7ZUuLa2lq1WhEzZzFuMK9AVVtrgVNHZV9oB26Qvt6ktrbWgq8p+0Ljg2ULTQUQ6sykaevivOV8RWiRQjMY4ZJB1dWpzqKVQDC2trYMc4fp4RbnQp2dnZ9//vnnn3++DW9sY2PDsv9cepnN5iaLZbiHDWlMJlOTf96GDRueeeYZQggAMAyjUqluk7gxhmECAwN79+7NPQwNDVWpWlzyoFKpbvNbOeniDiyFgjpVB8f2vpTsC+1YEcT7gmU/o+wLjQ8WL7QEDRwrJv00FnxJ0bFIoaXrINoTHNQy+cZyceT2eFkU5O/vf+PGjYaHWq32gQceaJzAw8PD3t5eq9VynUKtVuvv79/w2wsXLqSmpm7evJl7yDAM1wvs0aPHrYmbUKlUAwYMaGX8VqvVyrlO7+trPlZMJrq1d+Gc7Astrcg0/R6VWm3JcSHZFxofLF5o/f3Yn/KoWi7te7MsUmjHStgEP5kXVBO8rCceMGBAeXn5iRMnACAvL+/UqVPcdGBBQQE3aEkIGTNmzPbt2wGAZdmdO3fef//9DX++YcOGUaNGNY52999//w8//MD9/MMPPzROjFop3henCe+swgiXDbSXp4JmR5Qj3hfvxNIqqTo2XmHHzfPSI3Ryclq8ePG4cePGjRv3+++/z5o1i4tqO3bs+Oyzz06dOgUA8+fPHzp06PXr1/Py8kwm06RJk7i/NZlMmzZt+uSTTxq/4GuvvZaQkFBeXm4wGC5duvS///2Pj2zLW4KGzEnDpQJ3kKanvb2JLe43k6MwV8JSeq2KdsKVUC2jAGmF9KtByqoDfO2XnDVr1qBBg06ePDl16tT4+HjuyfHjx/ft25f7OSoqKisra8+ePS4uLqNGjbKzs+Oer6ur27hx4+DBgxu/WmhoaHZ29u7du21tbUeOHOns7MxTtmUs1ptkldJqEzgqZY9sW6ToqNKuhRWlny+ToqOPhuK/uEXZpdTLjmgchM6HdfHYKEZHR0dHRzd+xt/fv/GAp5+f39SpU5v8lZOT07Bhw259NU9Pz8cff5yPfCqEgw308CDpRXSgH7YCLUrVsy90V9a1sKJw+4geDRU6HyKmuK30AID3I1SUBNxNeFsU4FghjddgpZCtBF+SgufL3FaqXomDIljnFSTel+ApU7eRVUq97YmPvdD5QLyJ9SHZpbTaJHQ+RCxVjz1CJGv9NSRFj+tlWpSqownKuxZWFHsV9PQk6UV4Odi80jrIr6Y9PBRXCzAQKkgHJ2LHkEsV2Ao0D8/aVgIcHb2NFD3t40NUyqsEGAiVJV5DUnCasAUp2CNUADx39zYUuIOQg4FQWXCasCXFdVBQQ+9R3qCQ0vTXkCM6FutAs1J0Cl0spsTPrGS4cLQlqTqFDgopTYAjcVDhBEEzzBTSi2hfHyXWAQyEyhLtRS5W0AoL3IVCblL1rAIXyylTgganCZtxuph2ciYedkLnQwgYCJVFzUC0FzleiK1AUyk6Gu+L1UER8NzdZilz4wQHa77ixOOquVuYWEgvon0UOSikQAkanClvhjK30nMwECpOvIak4m7Cm50uoUFKHRRSoCgvcslAy+uFzofIKPmgXQyEipPgyxzVU1w215iSB4UUSM1ADE4Q3ExXAxX1tKu7QmsBBkLF8XUATzuSW46twD9SFXwtrEy4obaJIzq2n28rbuUuUxgIlQgP12giBc+UUZh4X5KqwwmCf6QqdQchR7mfXMnicbFAI/nVYKinXdwwECpIfw2TVogTBP9IUfbsAAZCJcJt9Y2l6Nh4jXIHhZTJ2x687EhOGdYCAIA6M5wupnHeyq0EGAiVqKcH0VbTolqh8yEOqXrcQahEuK2+QUYx7epOnNVC50M4WP+VSEUgzoccw1VzAKDUW3Ij3FbfAG9AhoFQoeJ9cTchAECdGc6U0FgFDwopFk4QNMAbkGEgVKh4XwbHhQAgvYh2U/agkGL19CD5OEEAAMo+U4aDgVCh4jXkeCE1Kb5PiFvpFYsh0MeHpCl+guBqJTWxNMRF0bUAA6FCudtCJyeSWar0VgCvhZUMdxMCQIqO9lfwDkKO0j+/kuGqOQA4ioFQweI1OEGAE4QAGAiVLF7xiwX+NFBKIVjZg0JK1s+XnChS+gRBiuKXjAIGQiVL8FX6+TJ4sprCudtCoDM5U6LcWlBjgnPlNEbxq6YxECpXV3dSYaT51ULnQzh41jZKUPbp22mFNMKD2KuEzofQMBAqFwHo60OOKng3ocLPV0Tw13oZ5QZCXDXNwUCoaPEaRrHThFUmOFdGY7ywFVA0hW+rx0ERDgZCRYtX8P2YjhfSSC9ip/hBIYXr4kYMRqqtVmItoABH9Ww/DIQYCBWury85VUzrzELnQwi4WA4BAAHo50uOKrJTeKGcOqlJByesBRgIlc3JBrq4kZPFSmwFUvUsLhlFABDvyyhzmjAFx0X/hoFQ6ZQ5R0IBjuppXx9sBRDEK3XhKB6r1AADodIpc9Xc+XLqjINCCAAA+viQ04qcIMAlow0wECpdvIYcUV4gxAlC1MDJBrq6kwyFTRBUGOFPA+3libUAAAMhCnUhFGhepbJaATxfETWmwFOWjuppb2+ixggAABgIEQD0U969CfGW3KgxBZ67m6JjcYKwAQZCxN2tXkGtQIURrlbioBD6hwI31KbqcFDkHxgIkeLux5Sqo7HexAa/++hvIS6EELiqmAkClsKxQtrPF+vAX7AgEMR6k9xyWm0SOh/WgjsI0a36+ijocjCrlPo6EB97ofMhGjwGQrPZnJWVlZ+ff5s0eXl5ubm5lDb9/uXn52dmZlZWVnIPjUZjaSNGo5GvTCuSnQp6epDjhUppBVJ0NB6vhdHNFDVNiMfNN8FXc3DlypXu3btPnjw5KirqpZdeujWByWR69NFHExISHnrood69excVFXHPV1ZWjh8/vkePHlOmTAkMDOTi6O7du319fcP+tn//fp6yrVjK2Vb/96AQtgLoJglKmibEs7ab4CsQzp8/f8SIEadOncrOzv7hhx8OHjzYJMH27dszMzPPnTuXnZ0dHh6+fPly7vkZM2YAgFarPX369PXr1728vLjn+/btW/K34cOH85RtxVLOepmzpdTPgXjjoBC6WW9vkltGK5Ux2IRb6ZvgJRAajcbt27dPnz4dALy8vCZMmLBly5YmaTZv3jx58mQnJycAmD59+ubNmwGgoqLim2+++c9//mMymerq6hwdHW1tbRv+pLCwsL6+no8MowQNSdGxSoiEKTpsAlAz7FTQy5OcKJJ/JSiqBV0N7e6OteAfvATCgoKCurq6sLAw7mFoaOjVq1ebpLl69WpoaGhDAq1WW19ff+nSJTs7u3fffTcmJiYgIGD69Olm818HH504cSI6OtrNze3hhx8uLS1t6a1Zlr127Vr6365cucLD55OhAEfiZEMulsu/FcDzFVFLFDJBkKpn+/kSFVaCRmz4eNGqqioAsLOz4x46OjoaDIZb09jb/zU+5eDgQCmtrq4uKSmprKwMDAw8f/58aWlpnz59Nm7c+PTTT/fv31+v17u6uhYXFz/00ENz5879/PPPm31rk8n03Xff7dmzh3s4fPjwBQsWtJTPhsU4CADivGyT8ur9gu9w5KLUC+1Igd2LYfUGg1XbO6kXmiCsX2hRLsy3f9q8FCrhYafWFNqBazYxbmAw1FghP2Lg6OioUt3hvqO8BEKNRgMApaWl3t7eAFBcXOzn53drmoaOXUlJib29vZubG5fsiSeeAAAPD49x48YdPnz46aef9vT05FJ6eXm9+uqrs2fPbumtbW1tX3/9dW6isTVcXFzu7rPJ18AA9mSpzbMud75TrXQLragWSuqNsR2cGatfDku30ARk5UIbEkxnpZucXVwk3Vm6Y6GdKDO9GaWS+Ke0MF6GRj08PEJDQ1NTU7mHqampvXv3bpImJiamIUFKSkpMTAwhJCwszMPDo3GAdHNza/KH+fn57u7ufGRb4RIUcDOaFB3b14dYPwoiSQhwJC625IKsJwhMLJwswhuQNcVLjxAAZs6c+dprrzk4OJw9e/bIkSNffPEFAFy8eHHo0KGZmZmurq7PP/98XFxcQkKCRqNZsmTJqlWrAMDe3v7FF1+cNWvWe++9d+nSpe+//55bbrp27Vo3N7egoKCcnJyFCxe+9957PGVbyaI8yVUDLa8HN9s7J5aoVD2N1+AOQtQi7qy1Lm6yjROnSmiQC5FxHW8bHgOhWq1etWqVh4dHUlKSr68vADg7O48ePVqtVgNA165df/755w8//LCmpmb58uUTJ07k/nDx4sXr1q179913vb29//jjj6ioKAAIDg7etm2bTqfz9/ffuHHjmDFjeMq2ktkwEO1N0grpfR1k2wqk6Oj8aAyEqEXcPqInuwidD97gDciaRW491UXSZs6c2blz51bOERoMBpy5aezN42Y7FVkUc7tQId1CM7Lg+bXx+iS19S+HpVtoAhKk0NKL6JMHzJkT+Ooh8O2OhTYp2TyyI5nWGS8Hb4LFgf4RryGpelboXPDlVDENxUEhdFuRnuRqJS2tEzofvMEbkDULAyH6R7wvc1RP5bqvHk/TQHdkw0Bvb3JcptvqtdW0ykTD5TsD2mYYCNE/vO1B40Cyy+TZCuBWetQa3ClLQueCFyk6mqDBRdPNwECIbiLjO5Ti4WqoNeJ9mVSZVgG8FmwJBkJ0E7nejEZbTatNNMwVWwF0B/EaklZIzTKsBNwNyLAKNAMDIbpJgi+R5eXwERwUQq3jZQcaB5JdKrdaUGeGzBIai1vpm4OBEN2khwfR1dDCWqHzYWl4AzbUerI8ZelEEe3uTpykujGEXxgI0U0YAn18SJrsWgG8JTdqvXg5jovgqunbwECImpLfbsI6M2SV0lhvbAVQq8jyfkw4KHIbGAhRU/G+jMwWjh4vpN3diSMOCqHWucedFNbKbYIgVc/GY4+wBRgIUVPxGpJeRE0y6hPioBC6K9wEwVEZjYtcMVACJMgZa0HzMBCiplzVEORMTpfIp1OI26fQ3ZLZbkLcRHt7GAhRMxI0stpWn6pjMRCiuxIvr4WjeC14exgIUTO4m9EInQvLuGygDCGBOCiE7kY/X5JeRI1yGRzFHuHtYSBEzZDTqrkUHe2PTQC6S65qCHGRyQRBlQnOl9NoL6wFLcJAiJrR2Y1UGemNKjm0Aql6iovlUBskyOXc3TQ9jfIidiqh8yFiGAhRMwhAP1/mqCw6hXhLbtQ2sjl3F1dN3xEGQtQ8ebQCVSa4WEGjcSs9unuyOXcXF4vdEQZC1Dx53I/pqJ5GeRFb/JqjuxfuRqpN9LrEJwgoQFoh7eeLdeB2sHRQ8/r6kMxSWmsWOh/tk4rjoqitCEC8RvITBOfKqKua+DsKnQ9xw0CImudgA13dSEaRtFsBPFYKtYcM9hGl4GKxVsBAiFok9ZvRUICjehwUQm0ngwkCPGu7NbCNQC2S+s1ocsuohx3xcxA6H0iy+viQsxKfIMCt9K2BgRC1KEFDjugkfLRGCl4Lo/ZxsIFubiRdshMEZfVwrYpGeGAtuAMMhKhFQc7EhiFXDFJtBfB8RdR+kj53N1VH43yIDTbzd4IlhG5H0osFcFAItZ+kq0CqnsUq0BoYCNHtSHeasKwerlfRnjgohNonQUNSJDtBkKqn8bhYrBWwjNDtSPdmNKk62gcHhVC7BToTtTQnCFgKxwtpHx+8FrwzbCfQ7cR4kfPltNIodD7uHg4KIUuR6CaKzFIa4Ei87YXOhxRgIES3Y6eCXp7kuARXzaXocFAIWYZEz93FVdOthy0FugMpHj1spnCiiPbFVgBZgkR7hKk6PFOmtTAQojuI15BUvcQWC2SW0ABH4mkndD6QLPT2JhcqpDdBgHdfaj0MhOgO+muYVB2V1vUwbpxAFqRmINKTHCuUUiXQ10BxHe3mhrWgVTAQojvQOICrLTlfLqVWALfSI8tKkNo0Yaqe7edLGKwErYOBEN2Z5OZIcFAIWVa8L0mV1G5C3EF4V7Ck0J1Ja9WcvgZK6mhXHBRClpOgYVL0lJVMJcAlo3cHAyG6swRJ9QhT9Gw8Dgohi9I4gId0JgiMLJwspn0wELYaBkJ0Z5FeJK+SltYJnY/WScUdhIgHEjpl6WQxDXMhrmqh8yEd2F6gO1MRiPUmaRJZNYe35EZ8kNC5u7hq+m5hIEStkqCRxmIBIwunimkcnq+ILE1C92NKxWvBu8RXINTpdFOnTo2IiJgwYcKVK1duTcCy7LJly6Kiou69996ff/654fn6+vp33nknLi4uOjr69ddfb3h+7dq1vXv3jo+P/+6773jKM7qNeA0jifUyGUU03BUHhZDl9fIk16toiRQmCFJ1NAEnCO+GDU+vO2XKlKCgoB07dnzxxRdjx449c+YMITf9Yz766KPvvvtu8+bNV65cmTx5clpaWrdu3QDgX//6V15e3qpVq9zd3c+ePcsl3rx58+rVq7///vuKiopHHnkkODg4Pj6eYnytTwAAHFlJREFUp5yjZvXzJY/rqVn0oTAFN04gfqgIxPqQND0d1UnUX7AbVbSOpWGuos6k2PDSIzx37tyhQ4dWr14dHh6+bNkyvV5/8ODBJmk+/vjjxYsXR0REjB079uGHH16/fj0AnDlzZseOHdu3bx84cGBERMSkSZO4xB999NHcuXNjY2OHDBny73//+5NPPuEj2+g2vOwgwIlklYo9EqbiqnHEmwQpHDd4REcTcLHYXeKlvLKysrp27eri4gIAKpUqJiYmMzOzcYL6+vpz587FxcVxD+Pi4rgEx44di4uL+/bbb8eMGfPvf//78uXLXIKzZ8/GxsZyP8fGxjZ5NWQdkthWj1vpEX/ifSUwQYAThG3Ay9CoXq93d3dveOjh4aHT6RonKCoqopS6ublxD93d3fV6PQBcu3bt6NGjPXr0WLJkyc6dOwcMGJCbm2tvb19WVtaQ2MPDg0vcrLq6uiVLlqxatYp7OGbMmOXLl7eUuKqqqsmALbqNGDfVIS0z3kO8hXatCurMtr6ktrJS6KzcDL9pbSDCQuvlRI7p1eWGSpW48vWPqqqqQ1rbd6PNlZVi77lajaOjI8PcocvHSyB0c3OrqqpqeGgwGDw8PBon4MJkVVUV93xlZSX3g6urq729/erVq21sbGJjYzdv3rxv374HH3zQ0dGx4QUbEjfL1tb2pZdeeuKJJ7iHzs7Ozs7OLSWmlN7mt6iJxEC69pzZyclJtIV2Rs8O8BPj/xS/aW0gwkJzBujgZPrT6BTpKdJIWG2kuRXMgE52jnwt/5AnXoZGQ0JCLl++bDKZuIfnz58PCQlpnMDR0dHX1/f8+fMNCYKDgwEgNDTUxcXFxuav/6G7u3tlZSUABAcHX7hwoUniZhFCvLy8Qv/m6+tr2Y+mZPd4kOI6qq8VaRMAeNY24p/IJwhOljI9PAhGwbvFSyDs27evt7f3xo0bAeC3334rLi4eNWoUAKSlpX3wwQdcmilTpqxbt45lWb1e/+23306ZMgUARo8eXVtbu3fvXgBIT0/PyclJSEjgEn/00UdGo7GiomLDhg1cYmRlBKCvDzleLN5Ig+crIr7Fa0S9rT6tiME58jbgJRASQjZu3Lh06dLAwMCnnnrq66+/tre3B4CsrKytW7dyaRYsWFBVVeXn59e1a9cpU6YMHToUAOzs7L755punn346PDz8/vvv//zzz8PCwgBg1qxZnp6eAQEBQUFBgwYNmjhxIh/ZRncUasha/8WGU6dOCZ2RpmpqajZs+ubM7s09nKWwzwtJVk+7st+3fLlz506WFd0kXHJy8vfffBVuzBM6I9JDKG+3XGVZtqSkxN3dvWGo81alpaX29vYODg6NnzSbzWVlZV5eXk0Sl5eX29jYODk53eZNZ86c2blz5xkzZrQmhwaDgVvailrj9z17J85aXBE71ev05vVLX3lw7ANC5+gvLMtG9h9yvsNQk5ntVXgw49BesS2ywG9aG4iw0CoqKiIShlyLeNypPC/RrWzXd18JnaN/LHp3xbrf0suC7/VKW3/0123h4eFC50hKeBxLZhjG29v79mmaXfaiUqlujYIA0LBwFAni442bKx5eA50ii7sOWvfFYvEEwkuXLult/epHzgOA/I2ZV69evc0sMkJtlpKSUtb5Pjp0RiVA2uoBJpPpNlf5VrZpy/ayF5OBsSlxdNu8fef8uXOEzpGU4L5L1FrhQR3UeekAwFw9ERbUUejs/MPHx4fVXwJjDdRVsYWXm72KQqj9OnTooNaeAdYEFTpbKqIoCAAaXx+4fhYAHK+fCA0UUfWUBBH9I5HILZ4359SUZ86s+KzcJXDxr18KnZ1/uLu7D3/q1Z0rE73tmeVL3xLbeBqSjYiIiJkPDfl0df9is92Md9cInZ2bvL9y1bBJz3kai0cOHfTYo7iK4u7wOEcoCJwj5JvBYJh+wjHGi7zWSyzDCbVmCN9q+mWESrS7u/Cb1gZiLrTtV9j3TrPHx9uI5wv3zEFzoDN5pXOVaAtNzMTSliEJWRjN/CfTbDAKnY+/fZzNJvgS0UZBJD8PhTAmFn67JpZexNVK+lMeO7MHtudthAWH7lp3dzI0gPkkRxTLx2vNsOosOz8av8nIegjAwhhmYbpYbseyJIN98R6Vh53Q+ZAsbD5QWyyOYVaJo1P4UTbbX0N6YXcQWdeDwYyJhV9F0Cm8VIHdwfbCskNt0cWNDA1gPs4WuFNYY4LVZ9m3ovBrjKyN6xQuEkGn8J1T7Iwe2B1sF2xBUBstjmFWnxW4U/hRDnYHkWAeDGbMFH7JEzIUYnfQIrD4UBt1diPDOjAfCdcprDbByjPm+dgdRAIhAAujmUUZQnYKl51iZ/ZQudsKlwNZwEYEtd2iaGaNcJ3Cj3PYQf5MBHYHkXDGBzMMCNYpvFRBf85jZ2B3sN2wBFHbcZ3CD4XoFFaZYOUZM84OImERgLeE6xS+fZJ9GbuDloDtCGqXRdHMWiE6hR9ns4MDsDuIhDcuiGEAfs6z9uXgpQr6yzX2JewOWgIWImqXzm5kuNU7hVUmWJVpfjMSv71IeARgfjSzKJ21cqdw6Ul2Vk/sDloGNiWovRbFWHum8KNsNhG7g0g0xgUxNgz8dNV6l4MXK+iv19gX78EG3DKwHFF7hbuSER2YD7Ks1ApUmWB1pvlNnB1EYjI/ilmcYb1O4dIM7A5aErYmyAIWxjBrs8wVVukUfpjFDglgenpgdxCJyNggxoaBXVbpFF6soL/fwMWiloRFiSwg3JWM7GiNTmGVCdacNePJokiE5kcxS6zSKVyawb7cQ+Wq5v+dFAMbFGQZC6KZNWfNZfX8vssHWeyQAKa7O3YHkehwncIfee4UXiinf2B30NKwNJFlhLuSMZ2YD/nsFFaZYO1Z8wLsDiKxWhDNe6dw6Un25Z4qF+wOWhS2KchiFkQza7N47BR+kMUO7cB0w+4gEqsHAhlbPjuFF8rpnhvsS7hY1NKwQJHFhPHZKeS6g3iyKBK5BdGqBSdYnnqFS06ys7A7yANsVpAl8dcpXJfFDsPuIBK9+wOJow0vncIL5XTvDdw7yAssU2RJYa7k/kDLLx+tNMI6XCyKJGJhjGphuuU7hYszsDvIF2xZkIXNj2LWWbpTuC6LHd6B6eqG3UEkAWM6EUcb2GnRTuGFcrpPi91BvmCxIgsLcyUPBDLrLNcprDTCB1nmt7A7iKRjkaU7hYsz2NnYHeQNNi7I8t6KYj7IMpfWWebV1max93XE7iCSktGdiLMadvxpmcvB8+V0n5Z9AbuDvMGSRZYX5krGWqhTWGmED7PwZFEkPYuiVYsyLNMpXJzBvhKB3UEeYfuCeLEwhvko2wKdwjVn2RHYHUQSNKoTcVbDD+3uFOaU0X1a9vnu2FbzCAsX8SLImTwQyKzNMrfnRSqMsC7L/AZ2B5E0LYpWLW53p3DpSXYOdgd5hk0M4svCGObjbLY9ncK1Z9nRnbA7iKRqVCfiYQfb29EpzCmjSdgd5B+WL+JLkDMZG9T2TmGFET7A2UEkcW9FqZa0o1O4JIN9rZfKGbuDPMNWBvFoQXTbO4Vrz7JjApku2B1EUjayY9s7hdllNDmffa4bttK8wyJGPApyJuOCmDVn77pTyHUH34jE7yeSvPlRqsVt2lOI3UGrwYYG8Wt+NPNJDltyl53CNWfZ+7E7iGRhREfiaQ/fX7m7TmF2GT2A3UFrwVJG/ApyJuODmbV30yksr4cPs8zzsDuI5GJB9F3PFGJ30JqwrUG8mx91d53CNWfZB7A7iGTkvg7Eyx62tbpTmFVKDxXgYlHrwYJGvAt0Jg8Gt3amsLwePsrG7iCSmwXRqqWt7hQuyWDnRKgcbXjOE/obNjfIGt5qdadw9Vnz2CCmM3YHkbwM70C87GHr5Tt3CrNK6WEd+xx2B60IyxpZQ6AzeSiYWX2nTmF5PXySw+LeQSRLC6JVS0/euVO4OIN9rRd2B62Kx8I+c+bM7t273d3dJ02a5OLicmuCoqKiLVu21NTUjBs3rnPnztyTv//+e0VFBfezj4/P4MGDAUCr1R45cqThDwcMGODv789fzhEf3opioneYXu6h8rZvMc3qs+YHAplQF+wOIhka3oH42MOWy+yksBYv9bJKaYqObhyksmbGEF+X3r///vvgwYMrKip+++23hISE2traJglKSkp69+6dlpam1+vj4uLS09O552fPnv3JJ59s27Zt27Zt+/fv555MT09/4YUXtv2toKCAp2wj/gQ6k0dCbjdTiN1BJHsLolWLM1hzy53CRRnsnF4MdgetjK/yfvvtt999993nnnuOZdm4uLht27ZNnTq1cYL169f36NFj06ZNAODo6Lh8+fJt27Y1/G3//v2bvGDXrl23bt3KU26RdbzJdQp7qnya6xSuyjSPxe4gkrVhHYifA2xtoVOYVUpTdXQTdgetjper7+rq6iNHjowZMwYAGIYZPXr0nj17mqTZu3fv6NGjuZ/HjBnzxx9/NPwqKSlp48aNp0+fbpy+rKxsw4YNO3fuLC8v5yPPyAoCncnE0Ob3FGJ3ECnEwhjVohY6hYsy2NewOygEXopcq9UCgEaj4R76+fmlpqbemsbPz4/72d/fv6KiorKy0tnZOTg4+OLFixcuXJg9e/aTTz65atUqAFCr1X5+fmlpabm5uc8999zu3bujoqKafWuTybRz5868vDzuYa9evSZOnNhSPuvq6mxtbdv3WRWnnYU25x7o8xN5vovZ2+6mluD/TpMHOkGAbX2dhe5rLyr4TWsDuRbaAC/wsyf/O1f3WMhNz2eVQUoBs74f254qINdCaw9bW1tC7jDOxEsg5N6V0r9aOpZlVaqmnX1CSOMEAMAwDAD8+uuv3JOXLl3q2bPntGnTIiMjR44cOXLkSO75WbNmzZ079/fff2/p3R0dHT08PLifXVxcuJdtFsMwt/ktalY7Cy3QGR4JgQ9zYWn0P1/NsnpYfx6OjAGGkee4KH7T2kDGhTY/El48SiaGgE2jz/f2aZjTE5xt2/WRZVxovOIlEAYEBBBCCgoKgoKCAKCgoKCh89c4TcOal4KCAnd3d0dHx8YJwsLCwsLCcnNzIyMjGz8/YsSIHTt2tPTWNjY2991334wZM1qTT7VarVbjEUZ3p/2F9lYMjfrB9EovdcNM4boz5geDobOHbKdG8JvWBjIutPsCoUOm6YdrzOTwv4LW2VJ6rMj87RAbdfuaZBkXGq94uXZwcHC49957d+3aBQBms/mXX37h+nN1dXWZmZlc/2/kyJFcAgDYtWsXl4D7FScvL+/y5cvh4eFNnk9OTm7Ya4GkqJMTeTSUWZ3510xhWT18lsPiUTJIURbFqJacZE1/N2wL09m5kYwDzg4KhK+CX7hw4YQJEy5fvpyTk2NjY/PQQw8BwOXLl/+/vbsPiqpe4wD+w+VFlFlxF2JfeFnGS5hmihCUFq1dVgFRBoGk4UYCAhpKkTaN1+4dKp0JBi2aUqPu2IsGk6ESJK9KgBr4VmSQI4oMtju8xIILLOyye87948zsdQDf7uiew57v56+zzKPznWeW87Dn9zt7nnrqqYGBAVdX15SUlP3798fGxkql0m+//bauro4QcunSpbS0tGXLllEUVVJSkpSUFBgYSAhJSUkZHh729vb+448/zp07V1lZ+Yhig3XsWDJjyVFT9iKB+0yy57J5nWKGLzaLAp+skNrJZpHiDuoff5vx+wDd3EsfVtrsFRHu+99C3UPX3t5eVVXl5uYWHR3t7OxMCBkeHq6vr1+1apW9vT0hRKfTHTt2TK/XR0VFeXl5EUKMRuPp06evXLlib2+/dOnSoKAg5r+6efNmY2NjX1+fVCpVqVSWJcDJsrKy/Pz87vPS6NDQ0JR3+sNdPKymZZ41z6YNG711z9WLm6PtbXsQ4p32f7D5pp3uplMaTLXL+rJ+F//dy2HrwodwUcTmm/aIPMJByAoMwkftYTXtP0fK07ftdBJJhU6CrsZjtr3VDe+0/4PNN62vr8/nuTV2Lm6GfvXZY18GByy+97+5F5tv2iOChRlgx65du6i3To1uKdPJA4/eefcTgK3KLdg3ptyq33zMnPLlW//axXYcXsMgBHaYzWZi70gIoRxnjY3Z4s2DAHc1ZjASx9mEEOI0e8wm75+dPrBLCdjxz+zMfxdE0JL5btorcbHY/QS8sz0z7XjEOsO1qhmd53MP7GU7Dq9hEAI7NqVuWBuxUqPRLF68GHc+AQ8pFIqrF09fvnzZz2+3SCRiOw6vYRACa2QymUwmYzsFAGtmzZoVEhLCdgrAGiEAAPAbrwdhYGBgb28v2ymmGT8/PwMW9h+EXq+fP38+2ymmGY1Gg49KD6q1tXXVqlVsp5iWeD0IdTrd7V/eBvdjcHCQ7QjTDE3TaNqDoihKp9OxnWKaMZvNQ0NDbKeYlng9CAEAADAIAQCA12xt16i7u/vnn39uea7F3Y2PjyckJGDv/gOZOXNmZGQknnl2/yiKcnBwUKlUbAeZToxG49jYGJr2QPR6fW9vL5o2wYEDB+bNm3f3Glv7rtHx8fH6+nq2UwAAACeEhITc8/tXbW0QAgAAPBBc4AIAAF7DIAQAAF7DIAQAAF7DIAQAAF4T5OTksJ3BGrRa7fHjx69everj4zPlw9BNJlN1dfWZM2dEItGcOXOsn5CDRkZGysrKWlpa5HK5s7Pz5IJr167V1NRcuXJl7ty5eC42o7Ozs7S0tKenx9fX9y43mbS1tbW1tSkUCitG466Ojo7S0tK+vj5fX187O7spay5dulRZWalWqz08PGbOnGnlhBzU3t5eVlam1WoVCsWUTVOr1RUVFW1tbSKRCL+e90DzQEdHh4eHx/r16yMjI/39/bVa7YQCs9m8cuXKp59+OjU1VSQSnTp1ipWcnDIwMODv7x8eHp6QkCCRSK5fvz6h4NNPP5VKpfHx8evWrRMKhWVlZazk5JTq6mqRSJSamhoUFBQREUFR1JRlGo3G3d1dLBZbOR43/fjjj2KxeOPGjQEBAdHR0ZMLKIpKT0/38vJ65ZVX1qxZk5+fb/2QXFNSUiIWi9PS0hYtWvTyyy9PLqisrHR1dU1NTU1OTnZ1dcU57e54MQg3b96cnp7OHIeHh+fm5k4oqKioUCgUer2epul9+/YtX77c2hG5Jy8vT6VSMafyzZs3Z2RkTCjo6uoyGAzMcX5+fkBAgLUjck9ISMiBAwdomtbr9d7e3jU1NVOWxcTEbN++HYOQsWTJkoMHD9I0PTQ0JJVKGxoaJhQcOnRo3rx5AwMDLITjJIqi/P39i4uLaZoeHBwUi8UXLlyYUBMVFZWTk8Mc79ixIzY21toppxVerBGWlZXFxcUxx7GxseXl5RMKysvLV69ezVz9i4uLO3PmjFartXZKjikvL4+NjWUuucTFxU1umpeXl+Uis1QqxSMpent7m5ubY2NjCSHOzs6RkZGTm0YIKSoqcnJyWrt2rdUDctHNmzdbWlqYprm4uISHh09uWlFR0aZNm/r7+0+dOtXf389GTG5pb2/v6OiIjo4mhMyZM0elUk1umlgsHhkZYY71er2bm5u1U04rtvYVa5NRFNXd3S2Xy5mXcrlcrVZPqFGr1cHBwcyxu7u7o6OjWq3m+TOj1Wr17U3r7u42m80CgWBy5djYWF5e3saNG60bkHM0Go2Tk5PljCOXy1taWibU/PXXX++9995PP/109epVqwfkIo1GIxQKLStYcrn8xo0bE2quX78+MjJy5MgRT0/P+vr64uLisLAwqyflEI1G4+bmZlkonfKclpubm5iYGBERYTab7ezsDh8+bPWY04ntfyKkKIqiKMtiskAgMJlME2rMZvPt+xqmrOGb23siEAhomjabzVOWJSUlKRSKrKws6wbknPt5F23ZsuXtt9/28PCwbjTuYk7TlpdTNm1sbEwgEDQ1NZWUlLz//vuvv/66dTNyzv007cSJE2q1+qWXXkpISLhx40ZVVZV1M04ztv+J0N7e3t3dva+v74knniCE9PT0yGSyCTVSqdTyhF6dTjc6Ojq5hm9u70lPT4+bm9vk3bYURSUnJ9+6deuHH36Y8sMir0gkktHR0eHhYRcXF0JIT0+PVCq9vaCrq6u0tFQoFP7888/d3d0jIyMZGRk5OTkTynhFIpHodDqDweDk5ESmahohRCaThYaGMqd+pVK5ZcsWk8lkb2/75647kUgkWq3WcoWmp6fHcvHG4p133iksLFy9ejUhZO7cudu3b09MTGQh6zRh+58ICSErVqyorq5mjqurq5VKJXPc39/PfMpRKpXMvgZCSFVV1YIFC/A3u1KptDStqqrK0rTBwUGj0UgIoWn6tdde6+zsPHr0KHMW4zmZTPb4448zTaMoqra2dsWKFYQQk8nErGyJRKKvv/5apVKFhYUFBgY6ODiEhYXNnj2b5dysUigUPj4+NTU1hBCz2Xzy5ElL0yzr9C+++GJ7eztz3N7eLpFI+DwFCSH+/v4ikYh5usD4+HhdXR3TtPHx8YGBAaZGIBAwv6eEEIPBgL9T74HlzTpWceHCBaFQmJOT8+abb4rF4q6uLpqmmXfJ+fPnaZo2GAwLFixYv359fn6+h4fH4cOH2Y7Mvq6uLrFYnJ2dnZOTIxQKmUbRNO3p6fndd9/RNP3RRx/Z2dklJiamp6enp6dnZmaympcTvvrqK4lEkp+fHx8f/+STTxqNRpqmm5qaCCETbqVoaGjArlFGYWGhXC7fs2dPTExMQECAyWSiabq+vt7BwYEp0Gg0Mpls27ZtH374oaen52effcZqXk4oKCjw8fHZu3dvVFTUs88+y7y7KioqhEIhU7B79265XJ6fn5+XlyeRSPbs2cNqXq7jy9MnWltbjxw54ujomJiY6OPjQwihafqLL76IiYlhdjcMDg4ePHiwv79/5cqVoaGhbOflhK6urkOHDhmNxvj4+IULFzI/LCoqeuaZZ3x9fc+fP//LL79Yiu3t7VNSUlhKyiF1dXW1tbWPPfbYhg0bmG9m6OvrO378eFpa2u1l3d3d1dXVSUlJLMXkltra2rq6Og8Pj+TkZGbjTHd394kTJyzvKI1G88033xgMhrCwsGXLlrEalisqKysbGhrkcvmGDRuY6wp//vnnyZMnX331Vaagrq6usbHRzs5OqVQ+//zzrIblOr4MQgAAgCnxYo0QAADgTjAIAQCA1zAIAQCA1zAIAQCA1zAIAQCA1zAIAQCA1zAIAQCA1zAIAQCA1zAIAQCA1zAIAQCA1zAIAWxHRkaGt7d3Z2cn89JoNIaGhgYHB4+OjrKaC4DT8F2jALZjaGgoKCjI1dW1sbHR0dExOzt7//79Z8+eXbp0KdvRALgLgxDAply8eHH58uVZWVkvvPDCmjVrCgoKtm7dynYoAE7DIASwNR9//PEbb7zh4uISFhZWUlLCPNsdAO4EgxDA1ty6dcvb21un0/3222+LFi1iOw4A12GzDICt2bRp04wZM7y8vDIzM00mE9txALgOgxDAphQWFhYXF+/bt+/7779vamp699132U4EwHW4NApgO1pbW4ODg5OTkz/55BNCyAcffLBz587KykqVSsV2NADuwiAEsBEjIyPBwcECgaC5udnZ2ZkQQtP02rVrz5079+uvv0qlUrYDAnAUBiEAAPAa1ggBAIDXMAgBAIDXMAgBAIDXMAgBAIDXMAgBAIDXMAgBAIDXMAgBAIDX/gsBYcj2oQhvQgAAAABJRU5ErkJggg==", "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: 9%|ā–ˆā–Œ | ETA: 0:00:01\u001b[K\rDiagonalising Hamiltonian kblocks: 22%|ā–ˆā–ˆā–ˆā–Œ | ETA: 0:00:01\u001b[K\rDiagonalising Hamiltonian kblocks: 38%|ā–ˆā–ˆā–ˆā–ˆā–ˆā–ˆ | ETA: 0:00:01\u001b[K\rDiagonalising Hamiltonian kblocks: 53%|ā–ˆā–ˆā–ˆā–ˆā–ˆā–ˆā–ˆā–ˆā–Œ | ETA: 0:00:00\u001b[K\rDiagonalising Hamiltonian kblocks: 69%|ā–ˆā–ˆā–ˆā–ˆā–ˆā–ˆā–ˆā–ˆā–ˆā–ˆā–ˆ | ETA: 0:00:00\u001b[K\rDiagonalising Hamiltonian kblocks: 84%|ā–ˆā–ˆā–ˆā–ˆā–ˆā–ˆā–ˆā–ˆā–ˆā–ˆā–ˆā–ˆā–ˆā–Œ | ETA: 0:00:00\u001b[K\rDiagonalising Hamiltonian kblocks: 100%|ā–ˆā–ˆā–ˆā–ˆā–ˆā–ˆā–ˆā–ˆā–ˆā–ˆā–ˆā–ˆā–ˆā–ˆā–ˆā–ˆ| Time: 0:00:00\u001b[K\n" ] }, { "output_type": "execute_result", "data": { "text/plain": "Plot{Plots.GRBackend() n=126}", "image/png": "iVBORw0KGgoAAAANSUhEUgAAAlgAAAGQCAIAAAD9V4nPAAAABmJLR0QA/wD/AP+gvaeTAAAgAElEQVR4nOydd1xTVxvHn4QpiIBksEEQVx0MGRXcOAsuxI2jKloHVm2L1gEdVqytonWhdY8qrrrwtahVEZG664S6J6KAslfye/8IskzI4Gag+X74g9yce85DuLm/e855BgsAadGiRYsWLR8rbHUboEWLFi1atKgTrRBq0aJFi5aPGq0QatGiRYuWjxqtEGrRokWLlo8arRBq0aJFi5aPGq0QatGiRYuWjxqtEGrRokWLlo8arRBq0aJFi5aPGq0QatGiRYuWjxqtEGrRokWLlo8aXXUbUCsyMzPfvHnj5OSkbkPEk5JCLi7EVubDxo0b1LIlEVF2NuXmkrU1k50LhUIWi8Visd5/KzeXbt8mT08SCunIEQoMZHLc8tHZSv3slEByMnl6Kvc/rjwAAFDgM3/27Pnly0WWlvYvX17r3r2lvr6+MsxTF3l5dPkytW9fcWTVKpo0iejdJXr1KjVpQkZGdPgwffopWVgo156zZx//9Renb1+jwsL8hg1fNW/uoNzxKnH8+PFly4yzsnyePGEVFsZ7eHQTHX/7loRCKeempSXr6HiXlFBeXtr06W/nzGkq46DFxcU//bTf1jbg4cN/+/Y19PR0k9dsAGJvYpWpm1/Zd8TFxc2bN0/dVojn3Dnq3JmyspQ4xPffU0gI5ecTEY0aRatWMdx/cXGxQCAQ+9atWxQYSKmplJVFs2fTN98QszlrAeSL/rA6xYwZ+Oqrv3bu3Jmbm6tuW+RGIBAUFxcrcKJQKPjll6+Cgz89e/boB6aChYXUrx/98UfFkZEj6euvqbi44hJdupQiI4mIkpPJ1ZVOnFCWMfn59P331LPn6Z9/Hte2bcJnn40/c+aUsgYTR0RExOHDexMTWY8fF6enT83KItGPsTGZm0v5efRozL17bx8/poyMVVevrpd9UH19/bNn90ya1GLVqnFCodzXp1AoLC0tld4OdZmtW7cOGzZM3VaIIS8PLi7Ys0eJQ6xYgcaN8eIFAGzahBYtUFDA8BAFBQUlJSWS3t20CU2b4s0bZGbCzw8hISguZmxooVCYk5PDWHeqwt9/PIs10cBggZOTV25urrrNkY+SkpICxq+hukxpKQYOxIABKC0tO3LqFNhs7N8PVLpEX7+GlRUuXwaAEydgZ4ewMBQWMmmJQIDYWDg6IjgY9+5h5849I0d+uXOnMu8v4njx4gURn6gnURN7e3e5zp09+3s2247N9qlf31mVX23R453UZlohVApffIHRo5XY/9atsLPDgwcA8PQp+Hxcvcr8KDULIYDJk9GjB0pLkZeHzz5DQADy8pgZui4KYWlpKY/nSQQimJrOPnbsmLotkg+tEFZGIMCwYejevULSSkpgaoqBA8teVr5EN26Ep2eZXmZlYcgQtGyJf/9lxpLjx+HmBm9vnD3LTIcKc+cOuNySmTPnHT16VIHT3759e+3aNcatqhmtEKqN+HjY2SEzU1n9HzgAS0vcugUAQiF69sRPPyllIKlCWFKCzp0xdy4AlJZi7Fh4e+PVKwaGrotCCMDa2pUog0igp9fnzBklPJsoE60QliMUYuJEdO5cZZXF3x8NG0IgKG9TcYkKhejaFcuWVTTevBkcDqKjIRQqbsbt2wgOhosLYmNr1Q9T9OuHxYuRm5sr1ARrZEMrhOrhzRvY2+N//1NW/ydPgsfDxYtlL3/7DT4+FUs3zCJVCAG8fg1nZ+zcCQBCIb75Bi1a4PHj2g5dR4Vw5cq1+vr2Ojr2zs4zGjfGjRvqNkgetEJYzjffwNsb2dkVR3bsAJuNpKSKI9Uu0dRUcLlVrvz79/Hpp+jXT5FHw1evEBYGLhdRUQyvsirMuXNwdERBgVYINQ8NFMLhwxEWpqzOk5PB4+HUqbKX9+6Bx0NKirKGk0UIAVy9Cg6nQpujo+HoiCVL9gwZMnn16vWC8kdoeaijQujo6EX0H1Euh9NlyZIXPB7+/FPdNsmMVghFfPcdWrVCRkbFkbdvYWiIKVOqNHv/Ev3+e/TuXaVNSQkiImBlhbg4WUfPz0dUFDgchIYys7jCFO3aYds2AFoh1Dw0TQj374ezMxi/ga9atdHS0pXLdTU13XL4cNlBgQDt2+O33xgeqzIyCiGAvXvh6Ij09LKXU6bsY7GCiJLq15/8/fe/KjB0XRTC0tJSMzNv0R5hgwbz4uLiLl2CvT0iIjRiXUsqWiEE8NtvcHEp80Erx80N1tbVW75/iRYV4ZNPsHdv9ZYnT8rkQSMUlnnEBATg7l0F7VcSsbFo06ZsWVgrhBqHRglhejqsrXH+PMPdZmZmcjheREVEhQ0aeL5580Z0fOFCdOmi3Dus7EIIYPZs+PmhqAgAhgyZTHSeCETp7u69FBi6zgnh7dvo3RtGRv319BYR7TAzc8vMzATw7Bm8vDBkCGNuRMpDK4SbNsHevswHrZwVK6CjI0aZxF6i58/DxgZZWdUbv3mDoUNr8qA5cQLu7vDyQkKCwuYri+JiuLggPr7s5QcphHU7jlCjmDiRPv+cvL0Z7jYzM1MgcCDSJzIwNLTLysoiolu3KDqaNm4kaXGiquPHH8nMjL7+moiofXu3+vV3EGWyWJvc3NzVbZpyycigadOoY0fq3p3S0rYvWWI6btzDevUO6OubE5G1NZ05QwYG5OdHjx+r21Ytktm3j2bPpmPHyNGx4uCTJ/TllxQRQc7OMnXi7U0BATR3bvXjpqa0Ywd9/TV17UrLllFJSemZM2euX79ORCkpNGgQjRlD06bR+fPk58fQ38Mcq1eTiwv5+6vbDqWiAk1WHpozI1y/Hq6uZfMhZtm+Xair283Q8FsTk1leXr2EQmFxMTw8sGkT82NVxtKyNZEVkbWfXw8ZT8nORosWWLcOAoFg3ryo1q27eXrO9vHJV2C3v07MCIuLER0NDgdhYdUnASNGYP78Kkeio2FrW8XbQtP4mGeEx46Bzy+LBayMoyPatBF/iqRL9O1b2NoiMVH8Wffvw8enyMSkk5nZlIYN+7Vu/a1GecS8T3Y2rK2rTGTr0IxwzZrf9fTsly+X7ruoFUIGePIEfD5jYUPllJYiPBzOzrh2rXj//v0HDhwQLVR++y369GF4rGqsWbOGqA1RMdEbIhvZNSklBVwuzpwpeykUYtAgRUIqNV8IDx6EszP8/XHzpph3nz4Fh4OHD6scPHoUfD42blSJffLzEQphcXFxSUnJ2bPg8cRE6YWHQ19fosdKDZforl1o2VJifoljx04YGk4X7SUbGbV9+VLW3Qe1MHs2xo6tcqQOCaGurh3Rs61bpXvVa4WwtggE6NwZP//McLfZ2ejTB+3bV3igiLh4EVZW1TfzGWf06NFEI0XfVSLvhIRk2c/93/9gY4MnT8pe5uSgdWusWCGfAZoshJcuoUMHtGiBmqOKIyIwfHj1gykpaNYMoaGQee9VdXxsQjh//mIOx9XMrI2x8bK//67+7uXLYLOxYYPE02u+RPv0wcKF4t9KTk62sBhOBKICPr+NJovK06ewsKgeDVWHhJDNdiJCcLB0x3WtENaWpUvRrh3DkXx376JFC4SGVn+ozMtD06aIjWVyLLH8+eczIj7RT0QziLpaWCAmBrLHQURG5pqYDOTzPTw8ej5//vzhQ1hb48QJOQzQTCF89gyhoeDzER0t/T+enw97ezHzjMxMdOuGHj3EuFSol49KCF++fMnh+BEJiQQNGvhkVs1/IRCAy0XnzjX1UPMl+ugRuFz895/4d0eMmMrheHK5rTdv3qWI9api9OiydBmVqStC+PAhWKyRRKE//yzdB1crhLXi9m1wOEhNZbLPhARYWVXJUlHO1KkYMYLJscSSlAQ9PfTocatp0+aenp5nz76xsQGfDw+PihDGmpk7dyGbvYoILNbRgQNDASQkwNJS4n3hfTRNCPPyEBUFCwuEheHtW1nP2rwZPj5iPHtFi94uLmXpgTSEj0oI79+/37BhP9GaB4/X++nTp5XfHTwYRkZSHH2lXqK//opOnST6dRcUFJQqKREGQ1y7BktLMVd7nRDC27dRrx5atRLOmzf/33+lJ7bQCqHilJTAywtr1jDZZ0wM+HwcP17l4PXr1wMCRrdr9zmff7NynK8yOH8eenoIDgYqhU9kZKBHDzRvDgcH+PtLT5gyYsQ0ojNEIHrq5RUoOrhyJZo3l1VFNEcIK0d33b8v97menvjjD/HvbtsGHg8HD9beRmb4qITwyhUYGAw0MRlrYTG6Z88q95CTJ8FmS1n3hgyXaGkpPDywdWvtjVUPPXqI39HQfCG8dg2GhvDxAbRxhCogMhLdujEWyVc+S7hzp8rxgoICa2s3omSiJAsLtyJleKa+459/oKeHzz6rGLo8jlAoRFQU7OwwYwa4XISG4uVLif0kJZ3ncDx1dH5js7ssWVKRI3/CBPTtK9MSqxqFsKSk5MmTJ6I/PDkZvr5o21bx6K5z52BnJ3Fuce4cbGwQFaWorYzy8QjhpUuwssLevcJz586dP3++8m29oAD162PIEOmdyHKJXr0KS8vq2/x1gpMn4eQk3g1ew4VQtKDVpUvZS60QKpfLl6unFqwNGRnw9xe/b3Tnzh0eT7S1Dj5/aCqz67CVuH4dBgboVSn8/f2A+gMHwONhxQqEh8PCAhEREms//ffffxs3bgwPv+bpWbHTWVyMTp0wb550Y9QlhHfv3rW39+Dz+9rYePTvf8/WVr7NUbEEB2PBAonvPn0KT08MHYr8/FqNUns+EiG8eBGWlhJT33XqBAsLJp/Vpk/HmDFymqhuBAJ4eIhJkSNCk4Xw9Gno6aFfv4ojWiFUFikpKRcuXGvVCjt2MNNhampNnoSFhYUcjjvRaaK/bW3dZfmnKoBIBbt3r3JQbGaZO3fKrL11C8HBsLOrSSqEQnz2GSIjK46IknRLWi2sdKJ6hHDQoC+I/iYC0d+tWk1iRJzu3weHg+fPJTYoKEBICNzdkZz8ZPHi6J07d6ll9+hjEMLERHC5Epejt28Hm41//pGpKxkv0bw8NGpUfbNDw9m8Gd7eEte6NFYIDx+Gjg5GjapyUCuESmHs2JlcbmC9ekN5vEGMXA2iSN4avLSzsmBldadjxwmDBn2hpOngzZswMIC/f/XjklKsvX2Lvn3Rvj3S0nD+fNnioSQ/mpcvYWODc+cqjoiSdF+4UJNJ6hLCfv3GESURgSipX79xTHUbHo7PP6+pgVCIOXNe6Oi4slgb6tefMXjwRKaGlp0PXgjPngWfj7/+Ev9uVhYMDDBtmqy9yX6JxsXBxYX5utlKoqAADg417QVophCKyoOEhlY/rhVC5snOzuZyfUWrlBYWw65fv17LDmNiYGkpxRVz+HA5vpwKcPs2DA0rltQrU0OuUdGWob09kpPL3EmcnODvD7EfSVwcnJyquMns2wdHx5p2GdUlhBcu/Kuj42pqGmpp6Vr7/285ovQc5QU6xLJjxw49vV9EV5elZVumhpadD1sIz5wBn1+RMPN9WrdGo0ZydCjXJRocLNOOgCYQFYWgoJoaaKAQbtoENhtffy3mLa0QMk9eXp65eVsiARE4nD4ptSiAVFKCyZPRqlX1DL/VOHAATk7Ml7MoJyUFhobo2FH8u1KTbh86BD4f69cDQFERoqPL/GjS0pCcnBwRERX/7sYzYUL1/DLz5sHXV2JmKXUJ4Q8/oH//nIsXLzI++tq16NSppgbnzp1r2DCIqIQoxdbWj9nRZeEDFsLTp8Hj1RTJumgRdHTk8wqW6xJ98QI8Xh2oT5mZCR5PSlSPpgnh8uVgsyVuw2uFkHlevYKFxa8mJl5crt+YMTMU7uf1a3TqhN69pcQSvHoFa+uKdGWM8/AhjIzg5SWxgSzVJ1JSqsT+Z2QgPBwNGpyoV68jUayZ2YBVqzYAyMtDs2ZVtgZrzr6mFiFMSwOHI0ewo1wIBHB3x/79NbWJiFhsa+vh7Ny1YcN/3891omw+VCE8dgwcDk6elNjg4UPo6uKXX+TrVt5LdPVqeHvX1vFK2Xz5ZfWyi++jUUL4ww9gs7F4scQGWiFkmPx8+Pjghx+QnZ1dLQ+FLBQUFAwaNNHOzqNz5zGOjtnh4dK/EoMGiZ/sM8LDhzA2RtsaV+BkLMOUnY3+/eHnV5H4LTBwKtE5IhCltW1bFo1x6RJ4PDx6VHFiTg5atRIfq6QWIRw5Et9+q8T+RS7psqRXTkgAj6fqijwfpBAePQoOR8rug7093N3l7lneS1QggJ8fw2HHzCLy6kpLk9JMc4Rw7lyw2Vi7tqY2WiFkktJS9OtXPfmsXHz77U8GBr8QgcVa263bbKntd+5Es2bKcql//BjGxtK//LLXIyyPMhSVY1y4MNrQ8EcisFjbBwwIK2+2cCHat6+SnOzBA/HZ11QvhKLYMtmzxihGYKCsMw+RF1XN24rM8uEJYVwc+HyJFT/mzPlBT8+BzW6uq3tIgTwVClyid+6Ay0XVJDYaxODB+PFH6c00RAinTgWbLb0Ij1YImWTKFPTsWatEyd27f050lQhED/z8atyMflfjN1mOTNdy8PQp6teHq6v0lnIV5gVw+DD4fKxbh6KioqFDJ9naerRoMdzJKfPZs7IGAgG6dMGiRVXOEmVfq1b4VMVCKBSifXtV1IW4e1emh24R+/fDykp1G0sfmBCKrkZJhbKfP3/OZjciyiR6xWI5ZcivhIpdonPmYNAgeU9SBcnJsLFBbq70lpoghCEh0NGRGOlYGa0QMsZPP6FVK7yrDC83QiE2b4aZ2SEDg+5E+83M+m7cuL3mU4KCMGeOgsPVzMuXMDVF8+Yy7VXIK4QAUlPxyScIDUVi4j8//LD4xIkTS5eiZUu8fl3WQFSyqlqo1ooV1bOvqVgId+yAm5uK9m++/BKTJsnaePduWFvj9m1lGvSOD0kIRW5cNUQEHj16lM0OFPnostndT58+Le8Qil2ihYVo1kyD8uqV06VLmdebVNQuhEOGQFdXYhhMNbRCyAw7d8LOrqKokLxcuQIfH7Rvj+vXkZCQMGvWD8ePS6nCsGULWrRQStTRq1cwM0PTprLe8RUQQgDZ2Wjf/qSubgeinWZm/das2Th7NlxdK5Lm7NmD5s2rZx2rln1NlUKYnw8HB8h/J1SQzEz5qldu2gR7e7nTnCrAByOEu3eDx5MSqFpQUMRitSZaTfSbvr6tAte5wpfoqVNwcBCcOpWsvCxR8nLgAJo3l3XFS71C2KsXdHXl+LZquhBmZmbW4HIiFArT09OlrlcoWwhFXteKVdx98wZhYeByERMjRz7SZ8+kf4cV49UrmJujSRM5FngVE0IAw4ZVd5aZMQPt2lUsvISE4Isvqpwiyr5WXtVdlUL43XcYPFg1Q5WxfDm6dpWj/e+/w9GxeplfxlGLEJaWln777QJ3917Tp89nJI/url2wtsa1a1KadeiAevVed+s2pE+fIU8V2rVT+BItLi62sOhWr944DqfnzJnfKdADs5SW4pNPEBcna3vVC2FaWtpPP/109OjRXr1gYCDfxrnmCmFRUdGgQYMsLCw4HE5wcPD7V39iYiKHw7GysmrYsGHz5s0vSv67lSqEN2+KKQQhIwcPwt4eISESy1tLonfvKgnJas+rV68OHjz44MFbc3O4uMi3zamwEC5YsNTQ8CcisFh/DBs2FYBQiLFj0a1bmc9kTg5cXKovEImyr+3cCahQCEWl5GuO5mSckhL5bj0Ali5FkybKLcisFiFcvHiFkdFMopeGht/NmiWDq0aNbN0Ka2vpu6pffgldXfHJH2RH4Uv077//NjUNIwKRkMfzUOwrxiAxMeLzaUhCxUJ4+/ZtXV1rFutLoo46Or/Ku02guUK4evVqd3f3goKCwsJCT0/PlStXVmuQkZGRnp4OQCAQhIeHt5Xs4688IXz+HI6O2LJF7hNTU9G9O1xdqyQVk5Hff4era/VivLVh586dbLYdm92P6BMbmxvyfuMUFsKioqLBgyfa2nrw+UNHjSqb95eWYtAg9OtXJsaJibC0rH5nF2Vfu3hRdUI4YoSY0qMq4MgRNGsm3/86MhKtWlXstjKOWoSwd+/Pia4RgeixmVn/b77B1q24elV83YOa2bgRtrbS91PXrAGbjd27FbO3AoUv0aSkJAuLUUQgKuJy1VyhPjdXes6j905RqRAOGTKSaIXo49LRcZT3dM0VQl9f3zXvoml+//33Tz/9tIbGBw8edHZ2lvSukoQwOxtubtWdG6WSn4+ICPB4MpUvf5+nT+XbOpIFC4vWRFeIQHTM2VnuZCUKC2E5OTlo0aJiE76oCL17IySkbC9w3jz06FF93XjvXjg6YtGiVVOnTlVszUp2zp+HlRWys5U6iEQk1XurAdFuq/xRrDKheiE8eBAczlZd3eFEF+rXnxAauiI6GiEh8PBA/fpwckJAACIiEBuLGzfE7y8IhcJHjx5lZ2f//jusrXHzppQRExLAZleswNcGhYVQKBQGBY3jcjvUq+fauvUW9TpgRkTIXetbxUI4deo0onAiED3X1/+AhNDa2vrvd2kzzpw5Y2lp+X6boqKimJiYiIgINze3ffv2Sepq69at/fr1u/eOx0xURSouRo8eYpK31syhQ2jUCAEBCrrVCIXo2RM//aTIuTXA4bQhukQEoqMuLmoQQgApKeDxKh458/PRoQMmTwaAkhJ8+il++636KXx+V6JAFmu+rq7NQ6XtjAmF8PXF5s1K6l46t26Bx5N7hvfVV/DxUYp4q1IInzxBcDBcXPC//2HTpu1BQRPWrNkgqOTEVVyMGzcQG4uICAQEwMkJpqbw8EBICKKjER+P9HQUFhZ6evbg8fqYmHhwOHukpgR6+hSGhmVFp2tPLRctsrKycnOL2rXDkiXM2KMAL1+Cw8G9e/KdpWIhLCgoYLN7s1itdXRsVq6UOx+BjEKoSyonOzu7Xr16ot+NjY3fvn37fhsA9+/fz8zMzMvLy8vLq6G3c+fOde3aVfS7nZ1dXFxcLc2bPt2AiLVoUWFurkztnz9nRUYaJCezlywp8vcXEJGMJ1bm99/10tL0Jk7MV+DcGujU6Yc9ez5nsezZ7EurVm3OlbP3wsJCXV1dXd1aXSTW1rR8ue6AAfpnzhRYWICIdu5kBQTUCw8vnTeveO1adufO9Tw9Cz75RFh+yqtXd4keAKzSUv3FixdHRUXVxgBJxMbqFhbq9+vH8GcuO3Z2FBhoEBlJCxcWyX7W/Pk0fbpBjx7sP/8sNDICg/aUvoPBPt+nqIiWLdNftUpv6tSSNWuK9fWJqE9QUB8iys/Pr9zSwYEcHKhXr7KXWVmsGzfYt26xr15lb9/OvnWLXa/ewawsb4HgO6I8Y+NOlpY9avhXFheTm5tx48bYsIGZ/ziAagbLha6uLlC8fn1Jp05GrVoV+vgIGLBJTmbPNhg+HDxesVwfSH5+vlAoZLFYSrOrCmvX6rHZR27ffmFpaUJE8t7EhEKhgYGB9HbyK3RtcXJy+t///if6/fjx446ONc12z507V69evTwJFb4ZXxqNjETbtjJFlQIoLi5LMx0RIVPeLEk8eAAOh/m46Zs3oaODH398e+zYMUkfYM0wMiMUMXMmunWrWDF+9QqffFJWmX39erRsWSVcRFfXniiNCEQhS5TzwJyXBzs7nD2rjL7lID0dXC7k9aKv5nnEFCqYEf7vf3BxQf/+DHjACoVYtmynvv4PRCAqsLWVUq+jVStwOExGJTG1jX34MOzt5faqqz2iNDcKpNRR8YzQxER6+tMa0Nyl0X79+v3www+i3xcuXNinT58aGj979oyIJMVRMCuE27fDyUlK1o8dO/ba2LhbWrp+8cWqTz5BQEBtA7wEAnTqJHfCX6mUlMDCQkq5A6kwKIQlJejYsYpDbFoamjbFqlUAMHgwZs6seGv58lU6OrZstoOh4cQvv2Rk/OpERGDoUKX0LC+LF6NvX7nPKi3FkCHo149J1yqlCuGzZwgJgbMzDh9mrM/8/PxWrTpxucO5XJ9167bW0LJ/fxgYgImdkwoY9Of65hv06qXqfNx9+uDXXxU5UZVCGB4OQ8NapfTSXCE8duwYl8v9+++/T506xefzjx49Kjru5+d3/vx5ALt27dq9e/e///578uTJ7t27d+vWTVJXDArhyZOwssKdOzW1KSoq4nLbEOUSlbDZHX//XdEw+0osW4ZPP1XEuaZm/P1haqqI611lGBRCAGlpsLPDkSMVRx49gqMj1q9HVhYcHavUihPdZbKy4OHBjGtDZZ48AYej9LA8GSkqgouLrGkyKlNcjMBADBxYq9tEZZQkhKKFEwsLhIcznyZCIBBcv35d5GQuCVGBghqqTygGg0JYUoL27Zl3ERBLXl7eL7/8NmzY97a29xVbUVCZEBYVwcBAYn0lGdFcIQSwefNmPz8/Pz+/TZVypo4YMeLatWsA4uPjg4KCPDw8OnfuHBER8UZycjOmhPD6dVhaSlkoy8zEkiUZenr+7wrzjqshwFFG7t0Djyf3yphUNmwAmy0x17DsMCuEAJKSwOdX2ZxPTYW1NWJjcfo07OwqPEfK7zKvXqFlS5lyAcvOsGHMi2tt2LsXbdoo8jBUVIRevTBqFDOTCWUI4enTaNkSXbuqKEvc+xw8KL1AgWIwG+Hz4gVsbBSMWpYLP7++BgZLiGLNzNxfKxSLozIhHDUKpqa17USjhZApGBHCp0/h4CAxrqiwEAcPIiQE5uYICECLFgNMTL4yMor85JMOsny+NSAQoH17LF9emz7EICquNlt6fQvpMC6EAKKj0aZNlfxq166Bz0dcHL76Cv37lx2sfJd5+RLNm5dtKNaepCTY2sq6Dawy/P2xbp0iJ+bloWNHjB0rR/YiSTArhC9eICQEtrbq9Mu9dg26utVzGDEF46Guf/0FW1tZE7IrRmlpKY/nLXqUNzGJOKzQOrVqhDArC7q6DDzBaIVQJt6+RZs2WLq0+nGBAAkJCAsDjwdfX0RHl01WSktLDx48uJAy5XMAACAASURBVHv37trfMhYvhp8f8xsD1tbw8GCmK2UIIYCRI6uHLiUlgcfD8eNwdS0LOqx2l3n6FM7ODOykCoXw8sLWmraT1MOVK4oXgXr7Fl5emDattjYwJYQCAWJiYGGBsDCovKZkBRkZqF9fvpwpcqGMnA8REejcmfmNknKEQpibexPdIsrmcDqnpKQo0IlqhLB3b9jYMNCPVgilkJeXV1yMbt2qPzDeuIGICDRqhBYtEBGhrGTHt2+Dw6lee6j2BAfDyIixu4+ShDA3Fy1bVn/WO34cPB527waXi5QUMXeZx4/h5ITVq2s19JYt8PZmYPKkDMaMUXwe/+YNWrV6zOV24vE8Pv00IFuhMEOFhVAoFMbG7omMXHTr1q0LF+DpiU6dVFc9SiwCARwcYGenRA8UZQihQIDu3ZW1bl9SgrFj0br1dTe33i4uvps2/aFYPyoQwocPwWYz41elFUKJvHnzpk2bLjxee2Njj86db4mev54+RXQ03N1ha4uwMFy+zLy15ZSUwMuL+VrVe/eCzcaxY4x1qCQhBJCaCh6vepWcP/+ElRXmzkXbtigqEnOXefgQjo4KLiECyMuDvb2qy77LTloaGjZ8ERGx+sCBAwrcaHr1+pzoBBF0ddd9+60iDgYKC+HUqXNMTCYR7axXry2ff3vzZvU/avj6okED5ZZZVlIWwJcvYWuLd/FljJGbi1690LcvA7W+VSCEPj5o1oyZrrRCKJGvv/6BxVpLBKKr7u6BmzcjIAAWFggJQXy8Kr7DCxaga1eGB3r5Evr6DG+HKE8IARw4AAeH6uFTW7fCzg5ubsf09Zvo6toNGDCq2lmpqbC1xXYp9RzFM28ehg9X0FoV8PLlSxMTVxZrjYnJ5JCQqfKe7u3dh+gxEYhOGRt/OXgwlixBYqIcjpoKCOGrVzh2DGZmnkRCIhDtmD+f6Ugg+Zk2Dbq60tOt1RLlpcM9fRpWVkwWsk9Lg4cHPv+cGQdjZQvh5ctgsSRWVJYXrRBW8OwZjh5FVBSGDUPLlmCzZxMdIgLRCxar2cCB2L+f4djkGrh5EzweHj1iuFtnZzRpwnCfShVCAN98g65dq++IrF8PFsuB6CpRAYvVYet7G3rXr8PKCpLz7onn0SNYWDAcScYsf/zxh4HBLyJHBnNzD3kLQe/de8DCor2h4c9crvuxY5djYxEWBl9fGBujRQuEhCAmBjdu1LRUKIsQPn+Ow4fx/ffo1w/29jA1RadOsLbuQ5RIJDAxGf/nnwfks5tpVq4Emy1T7fJaotS88D/9BB8fZoJE792DiwvCwxnoSoSyhbBZM/j4MNbbxyKENjYtEhMTKx8sz1L47bfo2BHm5qhXD9bWcHaGnR1MTKCrO4WoOdEsIi8HBykJKZiisLBw5MgwOzuPhg0nrl7NsM/itGnQ12fe30zZQlhaiu7dMW9e9eNE9iI9IPpZbL6FK1fA58u3hTB0KCIiFDdVBZw9e9bcfDCRgOiekZGviQnatUNkJBITZX2Qv3Pnzq5du55UTXdbXIyLFyFKZt2iBUxN4euL8HAcPIjK0XcXLlywsHBu0MBxx44dlU9/9gwHD1ZJ+Onri7AwbN5cIauPHj1q166PnZ2H2qvrnT4NNltF/+jaCOGbN282b9584MABgYQHE6EQffvim29qYR8A4MIFWFkhJqa2/VRGqUJ47BjYbCadJz4WIWSxduro2P/22+NJk+DjAysr6OrC2BjGxtDVhaUlfH0RGoqoKMTGIiEB9+7h+vUbXG5rA4PPzc19/vhjj2pMnT9/kYHBQiKw2SunTp3DYM/Hj4PNRmwsg12WoWwhBPDyJezsqk/viGyIIon2ELl5eS0We+Lly+Dzq4Tn18C5c7Cz07iQifcJD19gbe3evHnHy5evFBQgPh7h4fDwgIkJ/P0RFSVfuRyxZGUhPr5M2Bo2hJUVgoMRHQ0iS6IVRNuIbNaseSxqwOfD3LyK8ql9868GRDm1VVZjWWEhzMnJcXRsq6+/uEGDLwMDR0lqlpmJRo2wf7/iFsbHw9IShw4p3oNYlCqENjbo3ZvJDj8WISQC0Rw9vV0ODujUCVOnYssWJCVJKWH6/Pnz3bt336k5kQyjBAaW11176OcXxFS3b9/CyEhZW18qEEIA58+Dz6/yDGhs7EQURjSTKIDHyw0LE7+gl5gIHg/vCplIRCCApye2bWPQZFWTloZt2zByJKys4OSECROwdy+ysmrVZ2kp7t/Hhg0YPRouLheJ3N/NwqdYW//x3Xc4dAjPnjH0ByifvDw0bAg3N9WNqLAQHj16zNh4tujT5vF8arhNJyeDz1fQcX3LFvD5Ssmmqzwh3LQJOjoMp139WISQxXrNYrmeZDyBEtNER+9nsXoTHTIzC1q3Tv6CvxJo04aZaBuxqEYIAaxYgdatK6Lsr1275uDgVb9+s8DAZVwu2rTBwIEQmzb877/B46Hq0nh1Nm+Gj49GT2Xk4vp1/PorevaEiQl8fDBvHhISUFKCt2/fDhgwzsHBMzT062r/tdevceECdu1CVBQmTEC3bmjcGAYGsLVFhw4YNQrz5hUQWRE9JMogampkdLVnT2zYoKzCh8qgRQtwubXNKSgXCgjh1auYORM83nVd3c+ISonSDQzcak71t3QpPD3l9mCIjoaTk5SEkQqjPCFs2BCjRjHc58cihHp6TUNDp6vbEClkZMDFBdOnn5gxI+LIkaNMdRsRAV1d5oMRy1GZEAIYPbrKvLb8LnP8OLhcdOgAb2/xs/y//hITiVFOTg5sbRnzQNMoCgtx/DjCw+HuDlNT2NvPYLPXEwkNDOYOGLBy5kz07w9XVzRoADMzuLsjKAhff41Vq3D0KFJSqt9bR40aS2RDZO3l1aWgoCybkmhHMCYG8nruqAaBQNC4cTs224HNbqGnl6BiTyjZhbA8NMvODuHhuHMHkZG/Wlm5NmrkM2PGKS4XUVE1OTEFBSEsTHarMHMmWrVi0um0GkoSwgULoK+vlGy0H4UQBgUxtsyoJIqL0aULZs1iuNt//gGbzXwwYmVUKYT5+XBzq4iXr3yX+fdf2NujXz/Y2YmP7zx8GHw+rlwR89bcuXIX4K6LpKejUaPe72pXJTVuPOXnn7F7Ny5dkmNWl5qaernq55uTg+3b0bcvTE3Rrx927FBnppj3WbBgAYsVQgSix7q6zioeXaoQvn0LUWgWh1NTaNbdu2jfHt26SfQkz85G06Yyre0XFmLwYHTurNwHF2UIoUAAIyMmXVsr9fxxCGGDBl0/+yxElfWx5GXCBPTvz3CGi6IimJszvKv8PqoUQgD//Qcut2xXo9pd5v59NG2Kvn3B5Yrf/N+/H1ZWuH69ysHHj8HhMB+popksWbLa1HQ40WFzc/+4OLmzKnTsGMRmN2Oz3eztPd9/980bbNmCzz6DqSmCghAbK36lWpUUFKBly2+JlhKBSKijYy9vD5s27bS19bCza7tjhyLBFpKEsLQU8fEICUHDhggIwObN0mPYRUnpLC0REyNeLK9dA5crJXF5VhY6dkRQEPOTqmooQwinTIGxsVLSAH0sQmhklNGw4dhz586p2xbx/Pwz3N2Z91f08wOHw1j9HUmoWAgBHDoEe3ukp4u5y2RkwNcXvXrBygorVog5d88eWFtXCaMODsb33yvZYk1i374D06fP//vv0/KemJWVxWY7EwmIwGJ9Nn/+8aQkpKbi9evq9+WsrLJZToMGZXf53Fw8ffq0desu5uatpkz5irE/RjKPH6NPH5Fz+HMWy5koisUKdHfvLlcnWVlZHI4nUR5RroWF+5s3ck91379Eb9xAeDj4fHh4IDpabqePmzfh5YXu3cU/uq1ejVatJD5/PH8OV1dIcitjFsaFMC8PenrM12QV8bEIoa5uAZs9adq002p/RH2fAwdga4snDFQtrMKyZdDRwbVrDHf7PqoXQgCzZ6NLF5SUiHnczstDQAD8/NCsGcLCxOQm3rwZdna4exd5eXmnTxfb2al/4qLJlJTgzBnMm4e2bV8RNRX5MbJYg9zcEr290bgxGjYUVRyDiwu8vdGrF0aMQFgYIiMRFYWJE+HtDWNj6OsHEcUQpbBYHVbXMhtsjRw/Dnd3sFhwcICogNvdu3cnT568cuXKmk8UCvHoEeLjsXIlpk5F9+6wsUllsYaI/mQ2e4Cu7gMrK7i64rPPMGYM5s3Db79h3z4kJuLBAzFzrLVr1+vp2bPZ9i1bdn78GFFRaNIEzZohIqJWe/YlJYiKApcrfmoYEoLx48WcdesWHBxUFybLuBAOGgQul8H+qiCjELIAUJ1l27ZtkycfcHIqatlyf1yczogRNGsWWVmp2ywiIrp6lbp3p8OHycuLyW7v3KGWLWnBAgoPZ7JbsRQWFurq6urq6ip9pEoIhdSy5eqUlIVErMaN7VJSzlZ+VyCgyZPpwgWysCA2m3btIlPTKqdv3EhTpnxhaHg1Ozt/3LjZq1cPUaXxdYL79+n4cTp+nOLjydqaAgPJ359mzOh282YOkam5+dP09OtsNlvUWCikzMyyn4wMMb+/fElXrzYjukNERLv4/L937VrTsSPDNi9dSosW0atX9OmntGIFubqWHV+wIHrduh18Pm/79ujGjRuLDmZnU2oqpabSnTtlv6SkkJkZNW1KTZpQkybUrBm5uAgHDOhy924PFkvQtOnfFy8eT09npaXR8+f08iU9e0bp6fT0KaWn07Nn9PIl1atHVlZlP5aWtHSpg0BwjsiaaKSJyagxY/xHjCBPT2b+2Js3afRo4nBo3Tqyta04npdHXl70zTc0alTFwfPnKSiIFi+mYcOYGV0qeXl5RkZGLBaLkd7S08nKirZto6FDGemvOkKhUCAQ6OnpSWmnLCFWCVu3bu3Tp4/o8eTePYSFwdwcISHK8huWnRcv4OCAPxRM7y4RgQB8Pnx9Ge5WEmqZEQLQ1bUX+X2wWIOXLFnyfoOoKDg5YfhwtGpVvdB8UlKSsfFIIhAV8PltNHnzWJW8fo3YWISGwtERlpYIDkZMTPUwwVOnTu2TN3MdAMDa2o1oOdFNFqujldUhHR3o6aF5c0yfXtucdjk5mD4d9etDTw/9+1dfbExISDA3H0hUTHTF2tp/wgR07gxraxgZwdUVgwdj3jxs344LF8Rn387Pz9+yZeu2bdtlybCakYEbN3D8OLZsweLFoiyAZcmPvvhicq3+SHGUlGDBAvB42LChyvEbN8Dl4sYNlJaWAvjzT3C5OMqYH7pMMDsj7NIFjo5MdSaGj2VptFqu0fR0RESAy0VAAAMl2hUjPx8+PligSAEAKXz2GerXV91yn7qEkM12eJfEecE0CXX2Nm6EjQ1mzoSNTZXYiRMnTpiZTRc5UPB47qXKq+2mYQiFwl279kydOic+vqzMeUkJLl5EVBT8/WFqWpGbRtJNTOHqE2lpaW3b9rCwaB0ePld05NQpjBgBPh9EMDJCx47YulW+7avbt9GrF3R00KAB5s4Vf+769Rt0dFaIBMnY2GPlSsTH49EjVcSMNm/egcUaRLRIR8cmNTVVSaPcuAFPT/ToUeV5YtWqNAODThyOl41NRyurF5cuKWlwiTAohCkpYLOh1CDwj1QIReTmIiYGLi7w9cXBgyoNphYKMXQoBg9mctAtW7aamrYwMvJgsy+osoqQuoSwVasuLFZfou+I3P74Q2I81MGD4HLx/ffg8VCeILOoqKhNmy4WFhM4nIDp0yNUY7Am8PPPKxo0GEl01NS0x5QpR4ODYW4ODw+EhyM+XqaIbGYr1IvIyMCiRfDxgZER2GzY2mLECDExMK8qzfWOHUObNmCx4OQkJcHYvXuP9PTciHYaG38zdOgkZi2Xys8//zx27Ni7ygvjBfBu15DHq9g1HDNmOot1gAhEhwYMqHU5ZvlhUAjd3NCmDSM9SeSjFkIRAgEOHoSXF1q1QkyM0r2KRcybB19fJmtZZGRksNl2RDeI/mWxnN8qtcZaVdQlhABWrVo1ZcqUvXtf8Hg4dUpis+TkshKGdnZYuLDsYGlpaWJi4q1bt1Rjqobg4dGb6GV5HOGePXLnhVGGEFbm2jWEhsLWFiwWDAzg44NFi/DPPzf09e3Y7JZ6erZffHGRzwebjY4dIUvt9Jkz0anTf3PnLtyyZbvqp/5KrT5RjRs30LYtevbEkycICBhDdJkIRNd69hypGgMqw5QQJiaCxaoe8sQ4WiGsICEBAQGwskJERFmGxqKiohs3bihWyLsGYmPRqBFevmSyz6NHj7LZge/c23r99ddfTPZeI2oUwvK7THw8OByclhwRcPcuXFwwYwbc3TFmjEqTbGkOr17ByWkmixVDVGJkNGvJklUKdKJsISynoAArVpSViCKaQbSPCER/s1gjJk6UNWZfVM/y9Wsl2yoZVQohKk0Nv/76rIlJa13db01MWp85o4Yy00wJobMzOnasfTdSkFEI2Urx1NEw/Pzo0CE6epTu3ydnZxo37rmTU7suXX52cemUkJDI1CgXL9LUqXTgAPF4THVJROTn5wekEMURHWGzr7dv357J3jUef3/asYMGDaLkZPENnJ3pzBk6dYpataLMTOrShZ4/Lzl16tS///6rWkvVRlwcublRx44RQUFXnJz8xoxhhYWFqtuomjA0pMmT6exZys0lS8tkIj4REfHNzK6tXk3160vv4b//aPx42rOHLCyUbKvGoKtL4eF04gTt2pVSUNC0tNQcaPrvvynqtktB/vyTHjyg7dvVbUc5SldkZaJAhfoHD+DuHkm0mwhED7y8Ahix5MEDWFsrxX1r0yaw2Xft7PybN++QpFr/H02YEYo4dgx8PpKTJbbPyUHPnujXD7NmFerrdzQ1DbOw6Dd58reqsFV9vHmD0FA4OdW0eiw7KpsRVmbnzp1stgOLNYPNdl6+XEo4oAhRNr61a5VtmhRUPCMsp2vX4UR3iED0X5cuQ1RvACMzQj4fAwYwYo4UtDNC8Tg6UteuLDZbQEREpSkp7CNHSCisVZ/Z2RQYSLNmUc+eTJhYiQcPaPx4+vpr58eP42/dOu3j48PwAHWE7t1p0ybq04cuXBDfoH59OnSILCzowIGzOjpt375dlpGxLzb2mEAgUK2lqiMujlq2JENDun6dGI/bUxmDBw++efN/CxZwrl37c+rUSbKcMmkStWhB48cr2zQNpWNHt3r1thJlsFhbGjZ0V7c5irBqFWVk0MaN6rajMirQZOWhwIwQQFpaWqNGnjzeEEtLtwULktu1g7MzoqIULPAmECAwEOPGKXKu1J4tLeHlxXzPMqI5M0IRcXHg82sqTisUYuzYJBarDxGI8lis1uPHCy9cUK6pqufNG3z+ORo1kl6LUS7UMiOUl5ozjakSdc0IS0pKvvrqO1fXHqNHf2dvXyw23aBSqf2MsEEDTJzIlDlS0DrL1ERJSUlqamr+u2y4Fy+WJckNDZU7GH/aNHTrppS0nz16oH59deb71zQhBLBvH3g81BA7derUKV1dH6K2RK1NTDwWLhQ6O6NFC0RFISNDiQarjL/+goMDQkPBtKdXHRDCK1fA40nJPa0y1CWElXn4sCzdoCojxGophHPmwNBQ6XmSy9EujdaErq6ui4tLvXr1RC89PGjLFrpxg6ysqH17Cgyk48dJltxz69fTsWMUG0uMpyFbt47i4+nYMZl8Bz4e+venVasoIIBu3BDfQEdHp0GD1kRnif4hKv36a/z3H8XE0P371LgxDRpEx4+r1mLmyM+nWbNozBhasYJiYsjERN0GqZasLAoKotWrqVkzdZuiMTg40LlzdPkyjRpFJSXqtkYGiovpl19o1izmb5i15CMVQrFYWVFkJD15QsHBNH06ubvT2rVUUCCx/enTNHcuHTpEZmYMW3LrFn3xBf3wA7Vrx3DPHwBBQbR8OfXoQTdvinm3Xbt2np45PF6/hg3bcziTBwxgZ2aSnx/FxNC9e+TvT9OnU/PmtGgRvX6tctNrQWIiubrS/ft0/ToFBKjbGpUjFNLw4RQcTAMGqNsUDcPcnI4fp4IC6t2bsrPVbY1kdu/eY2DgaGjYpKRk4rx5tXPKUAYqmJwqD4WXRmUhIQHBwbC0RHi4mAoSd+6Ax1NKcqCSElhYoEMH5nuWFw1cGi1n167qRZcq8/Tp0zdv3pSUICIC9vY4c6bKuxcvIjQU5uYIDpZYLlVzyM9HeDisrfHnn8odSJOXRiMj0amT6tbTZEETlkbLKS3FxIlo25bhIGaxKLY0amDgSHSXSEg0JjIyUhmGiUW7NFpb/PwoNpYSEqiggFq3pkGDKCmJUlJSAgJGd+gwpFu3f378kTp3Zn7crl1JKKT4eOZ7/pAYNIiWLKFu3ej2bTHv2tjYmJqa6upSZCStX09Dh1JkJJU7kHp4lC2W+vvTzJnUrFlZZQMRANLT06EZVVnOnyc3N7p/n/79l/r2Vbc1auLECVq7lrZv17j1NM1BR4dWr6aBA+nTTyk1Vd3WiKO0VEjkRMQi8kxJua9uc95D6YqsTJQ6I6xMVhYWL4aDg9DAwJ0oiei6sbFHlmJupjXyyy/Q0cGVK4x3rAiaPCMUsXEj7O3x339SmqWloXt3dO5cvd6CiMoTxG3bHjo6tuXzezs6tn1YrbCFaikoQHg4rKywV5Hy6YqgmTPCx49hZVVTaiF1oVEzwnI2bgSfD6XWKVdsRujs3IdoBNESNtv28vvZZpWG1muUeV68SDc17SXKdsbhTDp//jyz/V+7Bh0dLF7MbK+Ko/lCCGD9etjbS6+GKhQiOhpWVoiLE98gKwvLl8PMbDJRPBFYrPhhw6bIaXVtOXcuafDgSeHhP5w8+bZ5cwwaJHeV89qggUJYXAxfX/z6q7rtEIdmCiGAY8dgaYkjR5TVvwJCWFKCBg3g5rZx7Nix15WdXbQqWiFkHqFQ6OzszWIdJjplZeXObKrSggKYmqJHDwa7rC11QggB/P47HBxw7570lqdOwdYWYWESU5L26TOW6B8iEP1jbDx2yhQcPYp3UTbK5e7du1yuD1Eim71WX3/grl2qGLQyGiiEkyahb18N3cTVWCEE8M8/sLJSVvIdBYQwOBgmJurZ4tXuETIPi8U6dWrPmDEnhwyJPXFimwmjDuwdO5K+Ph0+zGCXHwtjx9KcOdSlCz14IKVlx4505Qrdv09+fnRf3D5FZORkPn+ChcU0Pn/Cxo1TbG1p4ULi8ykggFatkt6/YhQU0NmzNH/++YyMwUTthMLxZmYPBw1Sylh1iJ07KT6eNm8mhmqhf0R4elJCAi1eTJGR6jaF6MoV2rtX07d4Ndg0jcTW1nb9+l8Z7zYigi5fpps3Nfpa0WTGj6eCAurcmU6dIkfHmlpyOHTwIC1fTj4+9NtvNHhwlXfd3Nzu3Dl58+bNTz75zszMjIjCw+nNG4qPp7g4+v57atiQevemXr2ofXvS11fc4NRUSk6m5GQ6f55u36ZWrahxY1cjo6m5uQNYrJv29lzFu/4gSEmhadPo2DEyNVW3KXUTUTL6zz6jp09pzRp13lh696bOnSkwUG0GyIQKJqfKQ8VLo0oiKQlsNtasUbcd71FXlkbLWboULi54+hRFRUVSV28uXkTjxggJkS9f140bVWq+R0dXqR4OICcn56U4H/bsbCQkICoKAQHgcmFlhYAAREUhIaFi6XX37v1eXoEDB4Y+f/5cDpsYQnOWRnNy0Lw5Nm1Stx01oslLo+WUJ6NncHlfrqXRqVNhaKjO9FgyLo2yoBlu4oqxbdu2o0ePbtegYh5yk59Plpbk70/79qnblPcoLCzU1dXVVcfDJIC8vLz68qfVWbqUIiOn6+md1dEpio6OHDq0pgDs7GyaOJGuXaNdu6hlS/kGysigkyfp+HE6eJCMjCgggAID6ebNDQsWrGKxzLy8bPbv35SSwrp0iRIT6exZevKEWrcmDw/y8KCOHcnBQd6/TOmUlpaWlpYaGhqq2xAaNIi4XFq5Ut121IjCl6iKKS6mzz+ne/fo0CHicBjoMC8vz8jIiCXDgvW9e9SkCa1aRRMmMDCuYgiFQoFAoKenV3Mz9Qjh8+fPf/vtt4yMjB49egQFBVV7t7Cw8MiRI4mJiYWFhd7e3sOHD5d0L/4AhLB1a3r1ip49I7bmbdfWRSH8559/Ondelp+/nSjf0tL3xYsrUk/ZsoVmzqS5c2naNEVMFQgoOZni4iguDlevugEXiPR0dcfr649u1MjX25t8fMjHh1q0IB0dRfpXGRoihEuW0M6dlJBABgbqNUQKdUUIiQig776j3bvp6FGyt69tb7ILoYsL1a9PV6R/BZWIjEKohrtvfn5+u3btsrOzfX19Z86cuXbt2moNTp48uWzZMmtra1dX18WLF48ePVr1RqqG6dPpzh365x9NVME6Sk5Ojr6+FRERGZWU6MhShmnkSDpzhjZsoJEjKTdX7hF1dKhdO/rxR7p4ERYWbCIdIjI0rLd7d8mNG7R+PY0fT61aaboKagjnz9PPP9OuXZqugnULFosiI2nCBPL1patX6e3bty9evFD2oEuX0sOHFBen7HEYQqnrs2LZsGGDp6en6PeDBw86OztXW3GuvC919epVHR2dfAkr3HV6j/D4cbDZ2L5d3XZIps7tEQIoLCxs0aJ9/fphRkb9jY1nyZ6aoKAAYWFo3Ph1ixY9eTwPT89er1+/lnf0uXMXcThduNwQb+/e6vroFEbte4QvX8LOTmKgp6ZRJ/YIq7F3L0xM1pmZefJ43Xv0GCoQCBToRJY9wlevoKeHuXMVspJRNDd8IjExsfO71GSdO3e+d+/ey5cvKzeovBaXmZlpbGxs8ME9H755Q4GBNHw4DRumblM+OFgsFov1WF8/q2PHlr160dmzMp1laEjLlpG9fdStW5+np1+8fHnM7NmL5B36hx++uXhxQ3z8zKSkw2pZUq67CAQ0YgSNHUu9eqnblA+X/v2ho7PizZvE9PRjFy40SExMVNJAvXuTsDJvhQAAIABJREFUjQ398IOSumceNXxX09LSmjZtKvq9fv369erVe/HihaWl5fstCwoKvvzyyzlz5rAlLx1euHChfJeRx+P9+ivzsQ3KwNPTkMNhrVlTkJ+vblMko949wvz8/Br+75I4derU06dtc3J+JcKFC203buwfFFRvzZriHj1kKlWvo/OCqBkRCQTNjhw5dOhQUadOArlWNblcLpfLLaihaommItojFApVWhkAwLp1W+LiEouLvYkmzJhRqsnfiMoofInWktLS0oiIxX///U/nzl6RkV9J3f2qjFAo1NMj0ep9QYFBbm5uvvwft+iUGvYI9+3TuXTJ4MKFwvx89VeZEAqFsnxEarjHGRgYFBcXi34HUFJSInZ/vri4OCgoqFWrVl999VUNvdnY2Ax+FwtWv359tW/1y8LUqazHj1n37ws131o1CqFAIFDg8zE3N2ex0omIqFAgEPTsaXDoEPr2NViwAKNHS/cL++qrMVeuTMjOHlS//q5BgxYsWmQwZgwNGIARI+Dr+4FHdqvFWSYmZvP8+RdzcmazWOtmzVppbKyQw5I6UPgSrSU//rh07Vphfv621NTlDRqsjYiYKdfpkyaFrFzZBbDKz889csS/e3cdea9q0V8tSQgLC2nCBPa4cXB1rUWYLXOIVnGlt1Pm8qx4pk6dOnHiRNHvz549Y7FYb968qdamuLi4b9++QUFBNW+01Lk9wtTU1P37wWYrvaQOI9TFPUKhEJaWU+rX9zQza21qGvvPPwBw6xbs7fHLLzL18N9//23btu2/d5m8Hz1CdDTc3GBvj7AwXLqkgFF1A7XsEQYGfk50jQhEj9u166/i0WuDuvYIu3UbTnSHCET/WVoOiYnB9euQa7PvyZMnN2/efP1a2K4dRo6EDDtoVah5j7BXL5iby2ePUtHcPcL+/fsfOnQoJyeHiHbs2NGlSxdTU1MiSkpKun37NhEJBIKRI0cWFxdv3779g9louXLlir6+bbNmIf37u/buffGjLamjbH7/neztf0tLO/P69eWdO4MDA+niRWrenBISaN06mjVLeg+NGzcePnx448aNRS/t7WnaNLp8meLiyNycgoPpk08oMpLu3VPuH/KR0LWrt47OOqIn9eqt7trVR93m1AH69euiq7uQ6EL9+j9169b10iUaPJjMzcnPj2bNokOHKDNTSg+2trYtWrSwsGD99RdlZNBnnyniLC2WM2fo2DHav78OusGrQJOrIRQKBw8e3KRJkz59+vB4vOTkZNHxnj17zps3D0BsbCwRtWzZ0uMdj6tl73hHHZoRNmnSnuggEYjOWFl5qNscmahzM8IXL8DlVilideQI+HxcvAgAGRnw8cHEibV6XBUIkJCAsDDwePDwQHQ00tIU702jUMuMcPVqgY1NtI9Pv/nzF9UtP1t1zQg3bRLa2m4dOPCL33/fUj4ze/sW8fGIiEBAAExN4eSEkBDExODGDTEpy4uLi3Nzc0W/l5Zi/Hh4eSE9XVYDJM0IBQI0bIi+fRX7s5SFRlefEAqFly5diouLy8zMLD/45MmTV69eAcjJyblXFUl/SR0SQjs7b6IEIhD9a2HRWt3myESdE8LgYMyZU/3g3r2wssK//wJAbi66d0e/fqj9Db+0FPHxCAmBmRl8fRETg+xsPHjwICBg9Kef9j169K/aDqByVC+EDx6Aw8GNG6ockzHUIoQZGbC0RFJSTW1KS3HjBmJiEBKCFi3KcgFGRCA+Hnl52LBhB5fbmsv9NDg4VKRnQiEiItC8OWSsvylJCEeNgrGxxLou6qJWQvjq1as9e/bMnz9/0qRJYWFhP/74Y3x8fJ5cORlVQh0SwpEjjxK1YbEmsNmN1q1br25zZKJuCeGRI3BxEa9we/bAygqiOmhFRQgORpcuYKqIVm4uduwoexI3MfEjOk10j8tt9+jRI2YGUBUqFkKhEN26aVD1TXlRixCOHYupU+U75elTxMZi+nT4+MDISKin14aogAgWFmPOVSrg+9tvsLGRqSS4WCG8dg1sNnbskM82FaCgEP711199+vTReecwXjmGz9jYeNy4cSouq1gzdUUIExPBZmPOnMcrVqxITU1VtzmyUoeEMDsb9vY4flxig9hYWFuXTT5KSxEaCk9PhsvepqeXmph4i+o2s1jzW7c+PHYsfvwRO3YgKQniEnFrFioWwuXL8emnKC1V2YAMo3ohTEiArS3evlW8h7w8gbm5O5GQCLq6M9auja/87r59sLTE6dNSOhErhLa2aNdOccOUh9xC+ODBg27duuno6HTv3n3NmjXXrl0rLCwUvZWdnX327NlFixa5ubnp6OiMGzeufIlZvdQJIczIgJERBg5Utx3yU4eEcMoUjB0rpc2uXbC2xs2bov7LloMk7D4riLt7dx2dbUSnudy2f/2VFhuLqCiEhsLfH05OMDSEkxP8/REaiqgoxMbi4sWKxPxCoXDEiKkcThsbG7cTJ04xaZZsqFII790Dj4eUFNWMphRULIQlJWjdGnv21LafmTO/43B6N2wYamfXlc8vCglB5UonJ0+Cx0PNRaHfF8KICOjpISOjtrYpA7mFcOfOnRMnTpTkllLOxYsXAwICUjTjEq4TQmhvDycndRuhEHVFCJOTYW2NStvNEtm0CVZWuHWr7OWyZXB0xJ07ilr5Hq9fv/7yy3nDhk25evXq++++fYsrV7BvH375BZMmoVcvNG0KAwNYW8PXF926HTc0HE8EonRHR2/GbJIZlQmhQIAOHbBsmQqGUiIqFsIFC9CrFzNd3b59OykpqbS0NCcHERGwsEB4eMVOwfXrsLfHqlUST68mhE+fQkdHc5e4NdpZhik0Xwh790a9ehr6rCSVOiGEJSVwdZXyDFuZDRtgb4+7d8tebt0KPl+K94FSEQrx5AlOn8akSbt0dRcSgUjI57ur3hKVCeHixfD11aBQM8VQpRA+fAguF/fvK6Xzx48REgJbW8TElK1UP3iApk0RHi7G4xTvCWHLlmjeXCmGMYKMQviBROlpJj//TMeOUVISNWyoblM+XBYtIisrGjRI1vZjxhBAXbrQ33+TkxONGEFmZtSnD23bRt27K9NQCbBYZGtLtrbUunX3uLhur18XlpTcNDAILCr6MCsw3LlDixZRUlIdDDVTH5Mn09dfU6NGSunczo62bKELF2jmTIqJoV9/pU6d6Nw5Cgigzz+ndetqqm7/++90+zb9959SDFMpUqXy7Nmzw4cP938HEyLNGJo8I0xIAJtdt9d/NH9GmJoKCwtFnpTXrYODQ8WJp0+Dx0NsrNz9MEtOTs6ff/6ZlJQ8aBD694eKP3sVzAhLSuDtjTVrlDqIilDZjPCPP9Cypdz5XxTj4EE4O8PfH9evIzcXvXohMBDVwgXKZ4Q5OTAwwNdfq8IwhWFmabSkpKR58+YHDhy4+A6GzGMGjRXCtDQYGmLIEHXbUTs0XAiFQnTtqvijxtq1cHTEgwdlL69fh50d1q4FAPVWIwJQVITevTFypErXD1UghAsWoGtX8QtudQ7VCOHbt7C1RWKissepoLgY0dHgchEaimfPMHYsvL2r+FeXC2GHDuDzVWeYYjCTYq2goKBt27Z9+vQpT/KimnlqnUYoJE9PsrGhP/5QtykfNOvXU3Y2TZ6s4Onjx9OXX1LnzvToERFRy5Z0+jQtXPicw/FzcPg/e+cd19T1hvEnYS9lbxQcqCg4QHCLiBt30WoFRxW3aGuNq0LrzxYcFZzFWhXBhVotqKg4UFwVRVHcggwVRIbsleT9/RFEVGZIchn5fvj4CTfnnvNETu5zz7nveY9jhw59JbB5aWXIy+PYMcTHY9EipiSInseP4euLPXsaee5y0bJ8OZyc0KuX5FqUk4O7O54+hYYGrKxgZARHR/Tvj8TEz4qdOYNr13D6tOSEiZVqjFBNTU1DQyNBcKmQUjOGDkVmJqKimNbRqHn3DitXws+vTju/u7vD3R2DBuHtWwAwM4O19fr09JWpqdeePVu+YsV6UakVAiUlhITgv/+wejWDKkQGl4upU+HlhRYtmJbScIiMxL//4rffGGhaUxNeXrhxA0+eICAA3bqhVy9ER5e+y+Vi0iRMnoxGMzKqyghbt26tqam5d+9eU1NTzY9ITFkD5fffcekSwsPRrBnTUho1ixZh5kx07VrXehYvxpw5GDCg1AsLC7MBAwBEhqmpWXWWWSeaNUNoKE6ehJcXs0JEwP/+B01NTJvGtI6GA5eL2bOxcSM0NBjTYG6OoCAcOYLYWCgpwcEB33/vpavbsVmz2UCBvz9jwkROVVGjsdIE+7UkLAyrV2PLlsZzo1Q/CQ3FvXvYt080tf3wA/j80jjSlSvd/vtvVmHh8JKSMy9fbk9MZHgEo62NixfRvz8UFLBkCZNK6sL9+9ixA1FR0knRWuDjA21tTJrEtA6gRw9cv45jxzBr1sE9e/4EfIBHRUXd2OwnTEsTGbUIYebxeA1x322J8e4dRo3C5MnCP7WSUhPy8jB/PnbsgJKSyOpcuhSTJsHBAWZmPR48OHX4cI/4+NPz59v17o3ISJG1Ihx6ejh/Hlu24O+/GVYiHMXFcHXF5s0wNmZaSsMhIQHe3tixg2kdH2Gx4OwMExNvYD4wBlhF1Ki84DMj7NGjx/bt2wWviei77767detW2btHjhxRVlaWqLqGA58PGxu0aoWAAKalNHZWrsSAAXB0FHG1Hh6YMAFDhiA5uai4uLiwsNDdHbt2wckJBw6IuK3a0qIFzp/HmjU4epRhJULg6QkzM3z3HdM6GhQLF8LdHR/3xKwv2Nv3A4KBPOAuUMK0HFHy2dRocnJydna24DURHTx4cOTIkT16SHfLrB5HR2Rl4dEjpnU0diIjceQIYmLEUvkvvyA+PsLWdqm8/Ehl5d+OH/9j2LA+Fy5g9Gg8fIjff2dyZq9tW4SGYvBgqKpi2DDGZNSWqCjs3Yv795nW0aD45x88e1Yfb3q2bt0aHNwlMbEtiyXj7b2YaTmiRJrdQQT89huuXsXly9IAGfFSUoKZM7F5M7S1xdVEWtoeLndXfv7qtLTd3t5/A7C0xO3buHEDEyeC2ScDVlY4eRLTpuHqVSZl1JyiIri6YssW6OkxLaXhkJODxYuxc2c9zSuUkHA/N/cFj5f4008/Ma1FlEiNsK6cP4+ff8auXdIAGbGzcSOMjcUbPqCvr8lmvwIAxLHZpTHS2toIC4OSEgYMQEqKGFuvlh49cOgQnJ1x5w6TMmrIqlXo1AnOzkzraFD8/DMGDYKDA9M6mhhSI6wTyckYNQpTpmDGDKalNHZevMCmTdi2Tbyt/P47p2PHLTo63c3Mtv73H2f37tLjCgrYtw9jx8LWFnfvildD1Tg44K+/MGoUHj9mUka13LyJAwewZQvTOhoU0dE4fLgxrJZpcEiTbguPIINMu3ZoTOtp6idEmDcPq1eLK+9wGbq6ug8eXCIiFov18iXGjEFkJLZuhbw8WCxwOGjTBsOHw88PY8aIV0kVjBqF4mIMG1aaN7wekp+PqVOxfTt0dZmW0nDg8zF7Nry9oaPDtJSmx5cjwlWrVsnIyMjIyMjJyQGYPHmyzEemTJnChML6y4AByM7G9etM62gC7N2LrCwsXCih5lgsFoA2bXDzJt6/h4PDpxnR8eNx/jwWL4anp4TEVMg338DTEw4OXya+qidwOOjZE+PGMa2jQbFtGxQU4OrKtI4myWcjwgkTJqSlpTElpWGxciWuX0dkJFRVmZbS2ElLw6pVOHOmTtnUhENNDcePY/162Nri2DHY2gJA5864cQNjxuD5c+zZA0VFSasSMH06srMxaBCuXq1f0SiXL+PECTx4wLSOBkVyMtatQ3i4NOcAM3xmhBs2bGBKR0PBx2c7h7ORx2Pz+b/v3Tuh7im+pFTLwoWYNk0E2dSEQzAj2qkTRo6ElxemTwcAQ0NcvYrvv8fAgThxgrEJQHd3pKdjyBBcvsxkIq7y5OXBzQ27dkn34KwdCxdizhx06MC0jqaKNFimFhQXF//44+/FxTd5vHtE/7O2Fs9yNinlCA2tF4mnR4xARAQ2bMDs2SgpAQBFRQQGYtQo9OrF5PrRX3/FoEEYPhy5uYxpKM8PP8DeHsOHM62jQREaiuhorFjBtI4aUFhYuHbt2pCQEKaFiJhPRvjy5cv379/X5JwnT55kZmaKTVL95dixRD7fHNAHmgGdr1y5wrSiRk5+PhYsgJ8fVFSYlgKYm+PWLbx7B0dHpKYCHweLv/0GBwcweGVYvx6WlhgzBhkZ+U+ePCkuLmZKyYULCA2FdF6pVgg6+c6djM2x15w3b94oKZl6ez8cPZpjYtKoZsM+GWF0dHSrVq1++umnR5Xc3xLRlStXJk+e3KVLl6wshhPzS5KUFPj6onNnrFljBmQCm4G/gYsdpBMZYmblSvTrh0GDmNbxkWbNcOIEhg9H9+6fVvJNmIB//8XcufD2ZkYVi4WdOyEjc9fQsI+9/W9t2/ZMSkqSpICcnBx3958HDPhu0qR//voL6uqSbLzB4+mJXr1EnzJQHEyYMAGYApwGol+/rtGoqaHw6Rnh+PHjVVRUli1btnHjxvbt2/fo0cPc3FxTU5PL5WZkZERHR9+8efPt27dDhw69e/euqakpc5olBI+Hy5exaxfOncOQIdi4Ebt3yyQkrOXx1gMkJ6fZoYM0+ZxYuHjx4tix8wsKuDIyPyYlzWVazmcIRoGtWmHoUPzxR2mMX48e+O8/jB6Ne/fiPnzwyMz8sGbN/BEjhkpMlYwMiorWFxUFpKZ2ZLFCfv99x44dv0usdRcX9zNn7EpKvlNQ+EFVVQ/oLbGmGzoPH8LfHw8fMq2jOvLycOIEYmJ0AS4AgM+wIFHzWbDM0KFDhwwZEh4e7u/vf+nSpX0f97lhsVgdO3acMGHCzJkzO3bsyIBMyfL0Kfbtw759MDaGmxv27IGKCtzccOkShgzp1anTkMLCwqKiWU5Oyhcu1JcghcbEiBGzioqOAe15vJFhYc0nT57MtKIvcXZGu3YYOxbXr2P7dsjKwsgI4eEwMHDNzV0PGE+fPjkyskPLli0lJqks2pBI0qngoqIelpTsAVBUNPXSpRu9e0uNsEbw+ZgzB//7X/1dbcnn48YNBATgyBE0awY2+28WqyNRFJCjoODGtDqRQpXz4cOHZ8+excfH5+XlVVGMQQICAiZPniyq2rKyyN+fHB3J0JA4HHr+/NNbS5eShgaNGkXFxaVH+HxasIB69qTsbFG1X+8oKCgoKSmRfLtsthlAAAEb5s6dK3kBNSQtjRwdacgQysggIuJyubq6dgLlSkprTp8+LUkxkZF39PW76up+p65ubWiYdOtWrWsoKSkpKCio7VnZ2WRgMIPF2gE81tAYev369Vo33DDh8/k5OTnCnfvixYvFi38ePfoPG5scHk+0ukRDTAxxOGRgQNbW9N13pKtLbm704QMRka+v77//3pOVpeXLmVZZA3g8XnHZVbtyqjLC+o/QRvjy5cs9e/ZERkYKfr1zh9zcSEODnJwoKIi+uPh7eFCzZuTkRF/8f/L5NG8e9erVaL2QKSPU0bECPIB/2GzT27dvS15AzeFyicOhNm3o4UMiImvrIbKy+4ALsrLWo0a9e/dOomLy80uDZU6fJj092r27dqcLYYTPnpGFBU2fnrNo0RpHx8lHj56oXZMNGaGNMC0tzcCgG3CSxdpqZeUkcmF1ISmJfHyoSxdq0YI4HPr3X+rdm2xsqPy3MDc3l8/n79lDbDZducKc1pohNcJKuXs3SlvbRkZmm7r6kPHj97duTRYW5OVFFV62tmwhVVUaOJAKCyt4l8+nOXOoTx8S9tawXsOUEa5endWihWu3bo6hoaGSb10IDh0iXV06epQyMjJ++umXadOW3Lnz0MODdHXJz4/4fAYkPX9OFhbk5vbl3VsV1NYIg4NJT4/+/lsYeY0AoY3wzJkzKio/C2YOdHV71uQyLW4yM0snw7S0yM2NIiIoJ4c4HNLRIR8f+mLMKjBCInJyIhWV+n7pkxphpcyezQHOAwRkaGk53rlTacl9+0hZmeztK3ZBAXw+ublR376Um1tbIfUdRowwI4N0den5c+HnnRjh3j0yM6M5c16NHDm9V68x586FEdH9+2RjQ0OHUnw8A5Kys2nsWOrTh5KTa1S+5kbI55OXF7VoQf/9VyeFDRrhjPDmTerTJ1ZGph+QATwwM7MTh7YqeP78uaPjt507D/rnn+DCQgoOJhcX0tQsnQwrKiIiCg6mFi3IxYVSUyuoocwIeTzS16fu3SX7AWqJ1AgrJjeX7O19WawNALFYl4YPd62s5PHjpKREffpU5YICeDxycaFBg6j2T1jqNYwY4U8/0dy5dXoAwxTv35OKSl8gHHipo9MzISGBiEpKyMeHdHTIy4u4XElLEjiWiUmNHKuGRpiVRaNHU79+lJIiAoUNl9p20QcPyNmZWrQgHx/av/9Y+/b9evce8/jxY/EprJC2bXsD/wEpiop91dVjBw6kvXspK6v03efPafBg6tKFbtyotIYyIySiuDiSlaUVK8SvW1ikRlgBwcFkakqTJxc4OEzW17fu2nVwUlJShSXPnSNFRerVq6bexuXSlCk0eHCj8kLJG+Hbt6SlRUlJDdIIuVyujs6nYJlTp06VvRUbSwMHUteuFBXFgLCQENLToz17qilWEyOMjqZWrWo349pYqXkXjYkhZ2fS0yMvL4avD1wuV0OjtIsqKHju2xdS9lZeHgkm8318qrljK2+ERPT338RmU0SE+FTXCakRfsbLlzR8OLVrR2Fh1Re+fp0UFcnWtna9lsulyZNp6NDqR5ANBckb4ezZxOEQ1S0kj0G6dh0sI3MAiJCTs+nZMyU6+tNbfD75+5OeHnE4DPSQZ8+oQ4dqDKxaIzx0iLS0aN8+0ctriNSkiz5+TC4upRaYny8ZXZWSlkaLFpGc3HAZmV3AWX1967JoruBgatmSnJ1rNMr/wgiJaNgwUlWtpw8La2iEjT/XaHExvL1hZ4fu3fHgQfUZHO7fx8CBsLDAlSu1S3okI4P9+6GujrFjUVRUF8lNlFev8M8/WLaMaR114Pz5AwsWPJk48dCtW3/NmaM3eDBcXSFIXMhiwdUV9+8jNhaWlggPl6iwsvxwAwfi3btan87jYflycDg4dw5Tp4pBXx24f/9+RkYG0yq+5NUrzJ6NAQPQsSPi4sDhQEmJMTElJfD1Rfv2ABAXd3DFivSZMy9cuhSgq6sbG4sRI/DTT9i9G0FBQm5jcuoUVFUbRnKcSpGAJ4uPakeEFy5Q+/bk5EQJCTWq8MkTUlIiKyvhb9+Ki2n0aBo3rjHMHUl4RDhpEq1bV/q6gY4IvyAjgxYtKg29Kz/dFBxMJibk5ibphTdlQS4VrkmpbEQoWC45dCilp4tdYa3IyclRU2vDZvdks1v+8MMqCbdeWReNjyc3N9LVJQ+P0oV3zHLhAnXsSI6OFBNDRPT27VtX10XDhrmGh1/38iItLfLwKI2RqSFfjwiJ6OXLevqwsL5Pjaanp1f2fK7mVGGEr1+Tiwu1aUM1j8CPjydVVbKwqOskRlERjRpF33zz5XrEBockjfDBA9LT+2QMjcMIBTx5QkOGUJcudPXqp4OZmeTmRkZGdELiS++Cg0lXl/bu/fJ4hUYoiIblcBiI9KmWqVNXAByAgGw2u+2BA3TnDtU8+QePx9u3b9/+/fuFa/3rLpqYSIsWkZYWcTiUmSlcraLkxQtydqa2bSko6NPBTp0GsFjBQJSMjJ2jY1INRwjlqdAIieivv4jNpvqWUKH+GiGfz58/f76mpmbLli1tbW1TvwrR5fP5Y8aMMTY2BnDt2rUqqqrQCIuLycen9I6s5g9jUlJIXZ1atxbNKoiiInJyImfnhu2FkjTCESNo69ZPvzYmIxQgCNRycvpsKcXVq9SuHTk7VxynLj6ePqX27b98ZPi1EQYGlq6PZJyiIrpwgdasoREjyNycVFWJxSIWaz/wP4CAQsBSWZlkZAggNpsUFEhDg0xNqWtXGj6cZs6kX3+loCC6d+/T6EdLqwOLNY7FGq2raymEJD6fn/nR7t69Iw6n1AIFaYaYJTeXPDxIU/PLa2BmJldNrSyeyzM4OFioyis2QiIaPJhUVWtxLyIBam2ET58+DSp/5yA2zp49a2JikpaWRkTffvvtwoULvy6zZ8+eBw8eaGho1NYIw8OpY0dycqK4uFpIysggbW0yNRXl896iIhoxgiZOrI+30jVEYkZ47Rq1bPnZN7bxGSER5eeTl1fpLVqZ4+TnU9nS+6SkpA0bfA4fPsIVf6fJyqKRI6lfv095JMobYUkJcTjUtm3pfJr4ePfunbZ2JxmZFurq7eM+fmnT0+nwYVq4kPr1IxMTUlQkgOTkSEeHbGzI1ZV27qTnz+ndu3cKCiZstjOb3WnixFlldT55QidP0saNNHcujR5NdnbUujXp6JCSErHZpU4pL/8cGPExk5/DgAFPv/+eOBz64w86eJCuXKG4uC8Xkpdn6NBxgAHQUkam9bJlPC0tWrSopus1xQqfT0FBpZEviYmfjguSZzVvTqqqA1isf4EoXV074ebkqjBCHo/09KhHD+G0i4VaG+GePXssLCwEr/X09MSXM9DV1XXZsmWC1xEREZqampWV1NTUrLkRJieTiwuZmNT6BjY3l/T1ydBQ9FFP+fk0cCBNm1bVN6o+IzEjtLcnf//PjjRKIxSQlFTaUct/5Pv3qXPnZAWFLizWHlXVHyZOnCMBJYL8cC1akCDVYJkRpqaSgwONGCGJ+T07u6HAdoCAg4qKk3R0SFaWAFJUJBMT6tePFi6kw4crfTxZUFBw4MCBe/fu1bxFHo8ePqTNm5+wWNYAH+ADXS0s3rZuTQYGpK5OSkokK0ssFgHEYpGsLCkpUfPmpK9PrVtT587Urx8BhsB7gIDZLVt61ZMllbdvU48eZGtLZZn4iFWYAAAgAElEQVRmk5PJx4esrKhdu9LkWW/fvnV1dR8+fOq1a0Je4aswQiJ68oRkZGj1auHqFj21NsITJ07o6+vzeDwSsxH269dv165dgtfJyckAsiuJGaiJEbZq1e3QoaAff0xUVy/57ru02l48CwvJxOSzp1OiJS+PBgygGTMapBdKxghPn6Z27b6cQ27ERijg8mWytCQHh9IkpUQUEHBQRmaTYIyio2MjMSWHD5OWFq1b99DKyqFt2z7bt4ebmhKHI4kem55OyspDgbsAAXGysv3WrKFLl2oXuyE0vXqNYrPN2GxTe/tvKisTH08REXT4MG3eTCtW0MyZNHYs2dsTYAzwAALWsdnzLS3J2ZlWr6aDB+nuXQaSTKWk0IwZZGBAe/YQj0c8HoWFkbMzNWtGzs4UFibKPH9VGyERbdlCbDbdvCmyFutCDY3w0zZMNjY2Hz58GDx4cOvWrbOzs9evX69XUSytn59fHeNUc3JylD6GEquoqADIyspSU1MTrrb4+NWTJ/+irLzexMT9w4diokM5OdWfVVhYePv27Q4dLB0djfPyWNHReQDV5EQhOHyY9c03SjNm8LdsKSzbK6dBUFhYKCsrKysrW31RYSHCihXKHh7FBQXcz49Tfn4+EYmvaWaxtsaVKwgIkBs0SGHMmJLVq4sNDfWaNTuWmbkIiM3IUBowgDd9erGTE1dOTrxKhg9HSAi7d+/hRBuA5vPnL96+PdjFRT0vT4yNxsSwly1TvHFDRlFxNuAGLGCx/GbO7Ll0aQ6AoiJJLEA6ezYwJyeHxWKpqqrmVPLl19SEpiY6d/7yuIGBQl7eYKArEHD27GF9/bwnT9hPn7KDg9lPn7JfvGDLyaF9e36HDvz27fnt2/PMzKhlS37Z1//UqbO//eanra35xx8r2rRpU5dPUVKC3bvl16+Xd3YuiYwsfv+etXq1XGCgrJYWTZ9esmULV0WFAOTm1qWRz8jPz+fxeKzKr2XTpuHECaUhQ2Revsyt1Qo0ccDn8xUUFOSq/RaVd8WQkJC+ffsaGBiw2Ww1NTWNiqi7RQ8dOnTrx7iI+Ph4FotVVMkdYE1GhLKyhUpKq8+cOVNzAQ8fPpSVNWSzhwAWKioXJBAUnptLvXqla2sP0dW17tx5YHJ9eJhQAyQwIjx0iLp3r+B2tdGPCMtIT6dFi0ozekybtlBJqZWOjmV4+I2gIHJ0JH194nAoNla8GuLi4thsG8FglMWa+eeff4qvrYAAatWKWCzq0IEEuXeOHj06duwEoaM3mWLOnDm9e/d+9OjR12/xeBQXR6Gh9Mcf5OZG/fuTjg6pq5OdHU2fThzOq2bN+gApwO02bXrVtt2SkpJVq37v3t1p5crfTp0qadeOhg2j6GgSdBhB1uzymRxETrUjQiIqKSFtberdW4wyakidokbFOjXK4XCmTJkieH3kyJGOHTtWVrImRqiqGqWj0yuhNlHAtrZDgEMAAS9UVNrX/MS64Oa2ksU6CBCLdeq77xZIptE6Im4j5HKpfXu6cKGCt5qOEQqIiiIbm+eysr2AO2z2nn79xgmOP39OHA7p6pKjIwUFiWtxKo/Hk5ExAu4CsWx26+flt+IUESUltHo1aWqSnBwNG1a7WLb6SW27aHo6Xb9Ou3fTN9+ckZX9+WOcTk9V1WI9PWrdmqytyd6enJxo4kRyc6OlS8nTkzZuJD8/OnyYTp2i8HC6c4d++GGTnNwPwFsWi6OltWnz5tK1qmLtIeWpiRHSx4eFZSuDmaJORhgcHCyI6qyMFy9e7Ny5UzhlL1++VFNTCwgIuHHjhrm5eVk948ePP3TokOD12bNng4KCVFVV165dGxQUVFmHCwgIMDKyOX36bM1bT08ndfVpwFHBYwkVFXPhPkVtGTfODbgNEPDM3n6iZBqtI+I2Qj8/cnSs+K2Ga4TJycnPnj0T4sQdO3ayWBsF10dZWauzZz89aiosLL3fNzAgDodevRKh3lKCg4M1NDqqqrbbvHmLaGtOTKSxY0lWlpo1oyVLGk8yXqG7aGpqqr5+N+CsnNxf3bsPy86m5GR68YLu3KFLlyg4mA4dIj8/Wr+e1qyhH34gNzeaMIGGD6d+/cjamuTkJgJPBffxCgp92rYlDw+Jbm9SQyMkIl9fYrOpiu19JIB41xGeOnWqXbt2wp1LRJcuXRoxYkTfvn19fHzK/k9XrVp14ePoYMmSJc7lqGwusVZJt/PyaMoUkpGhZs1esdkt2exxbHbLbduEtPPacv36DR2d7oqK3rKyPSdPrkHC03qAWI2woIBMTCp9ot5AjfB///PR0emrozO+X78xXy+B4PMpI4Pi4ujePQoPp3//pf37aetWWreOli2jHj02ArZALHCCzbZxcCAVFbK2Jg6HwsJK15Y8fVq6S5zg9l+0fxzhdqivgvPnqXNnYrGoVSuSyMosiVKXLvr48eOZM5euXv37h9rnnmnXrg8wHYgGvjcy6iX53S5rboRENGAANW8uoeinCqmhEbJIqHiE06dP//jjj0+fPhXiXBESGBgYGhp64MCBqotxuVi2DNu2QU0NPj5wcUFubu6FCxfs7OwMDAwkIxVAQkLCf//917q1zbRprSZOxOrVEmtZSMQaLLNxI27dwrFjFb9LRHl5eaqqquJoWkzw+XwDA+vU1DuAjILC/D59JsrL98vKguAnOxvZ2dDQQPPmn/00a1b6Ii0tYts2z8JCdUBdQeGhouLtfv3QujVYLDx8iMhIdO+O3r3Rpw/s7HD2LHbtwuPHcHHBvHlo0QIAiouL3717Z2hoKCMjI4R+LpfL5XIVhQpvKCwsTEpKatu2LQA+Hxs34o8/8P49bG2xYwe6dhWiyvoOI100OhojR/755s1xPl9DQeHDunVjfvxxniQFAMjLy1NWVq4iWKY8XC709WFhgatXxa2rYvh8Po/HqzZYRowBgfUBPh8eHti4EXJy8PbGkiWlx1VVVceMGSNhMS1btmzZsiWAixcxYABkZLBihYQl1BdycvDHH7hwgWkdoqOgAAEByMgQPPiBjAyvd2+Wnd2XhlclfXm8vocO/aumphIYuKN1a5w7h9BQnDsHfX1Mnw5jY7x7h+XL8fw57Ozg6IjZs3HnDqyt0aULhgx5sHGjS3GxcbNm727eDJHkHd7mzVuWLt0A6CsqZjo5Rf37bzMZGYwbh61boa4uMRWNnIICeHtj5078+uvsoiJ2aOj1YcMcFyyYybSuapCVxdWrsLLC77/X78udcOPNOk6Nioqqp0a9vUlZmRQVackSSYqqESkp1L49eXszraNKxDc1+vPPNH16VQUa0NRoSkppahhHR5o+fYeOTg9d3VGOjhN4IlqIx+PRnTvk5UW9e5OaGjk50caNtHMnubmRhQXp6tLYsTRlCikqfgPcBwg4PmnSXCEaEnpqVF7eFHgHELBeWdlb2OCBBoYku2hZKr6yNEAMUqupUQGCh4WhoW/Cw8NF9b2oIbVeR9iY+PNPrFiBggLMm4eNG8Guf5tN6ekhLAz29lBVxTxJz20wTFoaduxAZCTTOurMvXv480/88w/GjUN4ODp0ADA3I2NiTk6OYOgvEthsWFvD2hocDtLScPkyLlxAcDA0NTFyJCwtUVCA8HAUFT0CBLOaikeOKGVmwt4e/fqhe3eIY247NxeBgTh2DPfuobhYH1ABADQbPDh8zhzRN9dk+fABHA5On8bWrRg7lmk1wrJoEdav3zps2B42W1tB4cXr11GamppMi/oc4Wy23o4I9+whPT2SlaUpU5h8QltDEhLIzIzq7R20mEaE7u7k7l5Nmfo8IuTxKDiYHB3JzIy8vJjZnIjL/XKYqK1tA3QCZgLtVFUXKyqSoSEZGJCSEnXtSsuWUVhYNRGb1Y4ICwooIICGDSM9PWKxSEGBOnSgJUvo229/ZLMtWKxpsrJG4lh6UT+RQBdlaruuqhFiREhEsrItgByAgF/mzJknDmEV0lRGhDweT/AiJATz5iE5GaNGYd8+NGvGrK4a0aIFzp/HgAGQlcXM+j7bLxoSExEYiEePmNYhFFlZ2L0b27bB2BiLF2PMGAgVlSICZGQ+DROTkxEaiogIVaAEuAuoGhjkzZkDIiQnIy4OUVHw9cWff6KwEC1bol8/TJyIvn2/3He6uLg4NTXV1NS0/MHCQhw7hoMHERWF1FTIy6NVK0yejIULYWZWVmrj0qWTrl+/Pm2ab7MG8cWr96SkYMECPHiA/fthb8+0GtFAgGBqTobPr3cZo4Q0QjabLdbMWzXnyJHHp093NzGJfPYMQ4YgOhr1bcxdNW3a4NIlODhAWRmTJzOtRvx4emLePCE3wmaQly+xZQsOHsSwYTh2DNbWTAsqh4EBZszA5cuWhw715vFGysn5tmmjlJCAhAQkJiIxEfn5aNmyNG4lLQ1BQdi3D0TQ0oKtLYYNg4sLjIxa5eYWAzJsNuXlJVZofvPno3XrCgRcvXpt4sQFfL7GwYNnw8P/ES7uVIoAIgQEYNkyTJuGAwegoMC0IBExY8Y3u3d3IbIgynBwOMm0nC+pxfIJQYo5oZOCioPAwEAXlynAhDZtJl+7NqbBXV7LePQIgwZh61aMH8+0lHKIfPnE8+fo0wfPnkFDo5qSxNzyCS6Xm5ycbGBgIPjgly/Dxwe3bmHWLMybB0NDySuqERkZGRMnznv6NLZ3b5v9+33l5eXL3srPR3w8kpKQmIikJCQk4OVLxMYiLQ0A+HwQpQHdgJeAHDAE2CYvb25ggE6d0Ldv9X8sD4/+KSlHAH0FBe9t23Rnzpwuzg9aXxBHF42NhZsbCgrw11/o2FGEFYuSWi2fKE9sbOzjx48PHhz+zz8ycXEwMhKHui+p0/KJb7/91tHRcebns3XHjh1bunTpu3fvhPgvEDPaWVnnNmwY078/+vSp/ntbD+nYEaGhGDoU8vIYOZJpNWJj5UosW1av/0CxsbEODhOLikzk5BIXLjx68GArLhfu7jh8GB8TxddTNDU1w8IOV/iWsjIsLGBh8eVxPh8pKUhIgL//Az+/ZoDAO3W0tP4ePdpbcP8TF1d903l5xUAzAEVFmr/+mpeUhOHD0b17fQxSq7dwudi0CevXY9kyLF3K2JS7WGndunXr1q1HjkT79rC1RVJSfeohXz82LC4ulpeXF6SxTkhIKMtn/fbtWwBx9SlLYEBAAIvlzmLpX7iQ6uVFTk6krk6tWpGbGwUF0fv3TOurJffukb4+hYQwreMjog2WiYwkI6Oa7l7NVLDMhAlzgcsAAZcNDOaePy/K/WvqM4ABMBaYCeil1zL+5++/A7W1+2houLdoYR0WlurhQdbWpK1Nzs7k51cvtqsVByLsovfukbU1DRsm0UxpQiNcsEx5srJIRYWcnESlqCqET7EmMLyHDx8SUWBgYNu2bQXHuVwui8W6ffu2aIXWhYCAAAsLi/JpigTRdD4+5OxMmpqlpujvT18k5c7Ozo6IiEipJ/tpluPmTdLWptrspSFGRGuEjo7k51fTwowYIY9HXbvOBG4CBNwcM2amhAUwS8+ePS0sLGrrggISEhKuXLki2DlLQEoK+fuTszNpaJQmiouIaFR3FXXpomlpaVu37jx8+Eh2dokgZ17NvxqMU3cjJKJ79ySUklt4I0xLSwNw48YNIvL19dXS0hIcFxhkrLi3hKkNVS+o53IpJob8/MjZmbS0qFUrcnEhPz+6ejXO0LCLpqa7jk730NDzkhRcE65fJ11dunSJaR0iNcIrV6ht21qkxpe8EV6+TF26kJXVQy2trtrabvr6XR6W7ZnbNBB5rtGP1VJEBHE4ZG1NOjrk7Ez+/pSZ+alAVlbWjRs3MjIyRN60WBG6i2ZlZZmYdJOV3aGouFxFZbKzM6WmilydGBGJERLRzp3EZtOVK3WvqSrqlHTb0NBw3LhxDx48sLS01NPTO3LkCBEtW7ZMS0ursr0DGaHmSbd5PIqOpi1baPx4UlL6GfgXIOBNt27DxC1SCCIiSFdX7F2kWkRohD170sedRWqEJI0wIYFcXKhFC/L3Jz6fcnJy7ty5U29XMYoPMRlheeLiSm9M1dXJ2po8PCgoKEZfv4uW1gI9vW63bv0n1tZFi9Bd9OzZs2pqKwXbjKir96jJZbpeISojJKIxY0hJSbzPsOpkhHv27GGz2QAGDRr0119/AdDQ0ACwYcMGUeusE7XafaKMVat+k5HZBxAQo6s7rj499PxEWBjp6dF/jF4ZRGWEJ06QpSXVKrOSZIwwJ4c8PEhLizw8Gs/2QEIjASMsIy+PTp+mefNIRWU+cAUg4MHAgd9JpnWRIFwX5fPp118fsdlDgWLgrbFxN3FoEysiNEIiatOGTE1FVVkF1GlB/fTp03v27JmYmDhgwAA5OTl1dfV79+716tVrxIgREorhESc//jjnn39GZWQckpFJnTDB39YW48bB0xMSTFNcPY6O8PfHqFE4dQo2NkyrqQN8Pjw94eVVnyLEAD4fgYFYvhyOjoiJgb4+04KaGMrKGD4cw4ejqEhxz54cIgA5V64ojBwJJyeMGAFjY6YlioHnzzFvHjIzLebPdzpxoqeKivLu3VuZFsUw//0HExNMmICgIEZ1iNGLxY9wI0IB79+/F6R/TUsjDoe0tIjD+ezRRX3gxAnS0aG7d4nP54vwLqyGiGRE6O9PffrU+iyxjggvXqTOnWnAALp3T0wtNEgkOSIs4/Xr12Zm3fX0RrRoYf3gwYvgYHJzI319srAoja+RbIrmWlCrLpqfX5qZ3ceHvtqnsoEh2hEhEV2/Tmw2+fqKsMpPiHdj3npCXYzwC5KSyM2NdHXJy4vKxb4xz/HjpKa2TVPTSkfHcuNGiaYlrbsRFhdTq1YUHl7rE8VkhM+fk7MztWnTCPeJrTuMGCER8fn8lJSU8psSCGK/BcswdHTIxYWCgigrS/LSqqLmXfTiRTI3JycnSkoStyhJIHIjJCJvb5KRIXGsSKihEdan6SpGMTaGnx+uXMHdu2jXDrt2gctlWhMAoH//dOBARkbU+/f3vLz2ZWZmMq2oFvj5oUMH9O/PtA7gwwcsX46ePWFhgZgYODszLUjKR1gslp6eHrvc1Lkgk6qnJ+7cwe3b6NMH+/fD0BB9+sDbG8+efTqXx+M9fPgwJSWFAd01IDkZrq5wc8OWLQgJaZzzvSJh2TIMHAgHB2RnMyNAaoSf0b49goJw7BiCgtCpE44eRY0z0ImLrKwsJSUDQAaQKS7Wz2aqp9QGIvLx+bNPH+cVK7zWrClmVgyfj/370aED3r7Fo0fw9Gw8+RubAqamcHNDSAhSU8HhIC4ODg5o3Rru7jh9Ot/KymHgwI1WVuN37drPtNLP4POxaxcsLWFoiJgYDBnCtKB6T2go1NXRsyczrUuNsAJsbXHhArZtg7c37OwY3kjdzMzM3JyrobGwWbMFPB576dKWHz4wqacm7N17YM2aO9evry8oKDh2zFvCrefm5h45cuTixYtEdPEiunaFvz/OncP+/Q0v2beUMpSVMXIk/PyQlIQDB6Cqivnzg588Gfz+vf/79xfXrNlaT6ZwANy/j169cOAArl6Fl9eXu3xIqRA2G3fuIDYWU6cy0Hq92EGifuLoiMhIHDuG+fOhowMvL/Tpw4AMFot15cqJy5cvs9lsO7v+K1aga1cEBqJ3bwbE1JCwsFs5ObMAMx5vwcWLEu3XeXl5nTsPTE4eJSt7UV39mKLiznXrpBOhjQo2Gz16oEcPdOkiO3VqUUEBAG5GBktDA507w8YG3bvDxgbm5pB8UuS8PKxdiz178L//YdYsBgQ0aPT0EByMYcNgb4/pks3cXs2IMC8vb/fu3enp6ZJRU99gseDsjEeP4OYGFxcMGgRf3zOmprbGxtY7d+6TmAw2mz1w4MABAwYoK7N9fbF5M8aPh6cnPm7FWO8wM+vNYm0Hniorbxo8WKK3D9evX3//3qGgYFVOzq7MzKjoaK7UBRsro0ePsrK6q6vrpKPTLzBwZXIyvLxgZoZz5zBuHNTV0acP3N2xfz8ePar4GQePx+Pz+aLSExKCjh3x9i0eP4abm9QFhWHwYKxciVmzcP++ZBuuOpYmISEBQFRUlEgCeESOCKNGqyU/n377rZjNtgI+AIXa2n3evHkjmaa/JimJ+venAQPo9WsxtiJc1OipU6StTYsX7xk06Lt16zYLF3cqXNTo8+c0Zsw9NnsMwAcyDAy6CNF0k4WpqNE68vr169zc3K+Pf/hAERHk40MuLtSqFamrU+/exOFQcHBpHvC1azfr6HTW0bHy9t4uRLvlu+jr1zR+PJmb08WLdfgkDQRxRI1+Qd++1KwZiSRsXDTLJ6RGWJ60tDQdncGC3EiKim6hoZESa/preDzy8SEDAwoOFlcTQhjhkSOkq0u3btW16doa4cOH5OJCmprE4dDChWt1dTubmNicOVPvEsnWZxqoEdaQ16/p5ElatYqGDCFNTTI2TlVQ6A3wAG7z5nZXr2YmJFB2dk1ri4+P79ZtkJ6eta/vDh8f0tEhDw8qLBTnB6g3SMAIS0pIV5e6dhVBVXXKLCOlQrS0tNq0kcvL8+Dx1JSVYyZNspowgbGUNGw23N3RsycmTcLRo/jzTygrMyCjPIJcLWFhsLKSXKP37+OPPxAWhtmzEReH5s0BrN6yZbXkFEhpCBgZwcgIo0eX/hoenjN6tF5RERtAQYHOnDm5OTnqGRkoLoamJjQ0Pv1b4a+2tkOzs38CrNzdF3Tp0vvaNStzcyY/XSNDVha3b6NtW8ycid27JdJidYJkDQwMqt3et+kQHv7PP/+cKCgonDDhfH6+/KZNsLTEd99h5UpmIhJtbXHvHubOha0tDh2CpSUDGgTs3Alvb1y+jLZtJdTitWvw9sbDh1iyBH5+9X3jXCn1iv79zbp04cfEzAV43boph4WVLvErLkZGBjIzkZHx6UVmJp4//+zX7OwiYAYAYJKt7XZzcz8GP0ujpGVLnDyJkSNhb48pU8TfnggGn8wh4anRCklKokWLSjO0MbiTjL8/aWuTj48o66z51KiXF5mZkQh36Kp6ajQiggYOJDMz8vFpKvNREqBxT41+DY/Hu3LlSkREhBATfSoqbYATQByb3fXkyZPikFdvkcDUaBkcDsnK0rNnwtcgzSwjIYyN4euLqCgUFKBtWyxfDkbW+bm64to1+Ptj3DhkZEi0aW9v+PsjIgKtWom3ISKEhMDODnPnwtUVz5/D3V26Ol6KkLDZ7H79+vXp04dV+/jOW7dOmJl5qaqO4HBGji6bb5Uiary80LUrunYNNjHp2b//aMFeueJAaoSioUUL+Prizh1kZqJdO3h6MpArqF073LwJc3N064arVyXRIhF++AFBQbh6FUZGYmyIz0dICGxssGYN5s9HdDRcXSErfcAthSE6deoUG3szOfn2b7/9wrSWRs7ixUH5+d6vX/tcvdq3c+fBYmpFaoSixNQUfn64dQvJyTA3h7c38vMlKkBBAV5e+OsvfPcdli9HSYkY2+LxMGsWrl9HWBi0tUVZM5/PX7NmzYwZbjExMSUl2L8fFhbw9oanJ6Ki4Opav3Z0kiJFivg4evQ44AbYAUtTUsSVaVl6RRE9Zmbw88OlS3j0CK1bw9sbBQUSFTBoEO7exYMH6NUr3dZ2nJGRzbBhU/JF6sk8HmbMQGwsLlyApqYIKwaAtm17r1uXfPRon86dx7Rs+W7XLmzfjmvXMHKkdJGyFClNi6FDHVisvcBb4Cifb21vj3fvRN+K1AjFhYUF9u9HWFjpdha//57ZvftwPT2brl0HvRPHX/JzdHVx6hT4/LWRkZPevr1z6VLPdet8RVV5cTGcnZGejtBQqKmJqtZSSkoQH59CtBuYx+fP6dfP59o1DBwo4lakSJHSIJg9e/aMGVYqKgPMzP44ftw7NRWGhhg+XMTPnmpkhK9fvz5z5sybN29E2XLToFMnBAXh+HH4+W28e3dKauqd6Gj3H3/8nwSaZrOhpJQE2AIoLrY9dSpRJJny8vMxciTk5XHihCizCaenIzAQEydCTw9EfCAe4LNY4fb2LUXWhhQpUhogu3dvyc19Fhd3c9y41o8f4+RJPHwIHR388ANElSCvYiOcOHGip6en4HV4eLi5ufmIESNat2598uRJ0TTbxOjeHd26pRGZAyBqd+HC+9OnxfsAT8Dcud9qaroBe5s3X9KixcS2bTFhAi5cEH5vqawsDB4MXV0EBkIki0vj4uDri0GD0KoV/P3RqxcePsT+/b/Ly9vLyJja2mLOnDkiaEaKFCmNhZEjkZSEdevg54dmzeDjI4pKv15RUVJSoqCgcO7cOcGvPXr0sLKyioiImDFjRosWLbhcrvBrOkRNfVhHWEOuXr2mrW2roLBJS6v3jz+ec3QkTU1ycaGwMBLrmpzbt2/7+m57+PAhEX34QH5+ZGlJ7dqRlxelpVVz7hfrCDMyyM6O5s2jcnuJCwOXSxERxOFQhw6kq1u6//gXiwbFtEO9lCpoausI60iT7aKSXEdYNTweLVlCcnKkqUlHj1ZWRthco4InWC9fvhS8ZrPZBw4cIKK3b98CiIuLq5N2kdKAjJCIXr16deDAgRcvXgh+TUwkHx/q1o1MTGjRIrp7V3JK7twhNzfS0CBnZwoLq7RYeSNMSSFLS+JwatpEbm5uZGRkZmZm2ZGMDAoKIhcXUlcnCwvicCgiotKbgCZ7lWEQqRHWiibbReuPEQrIy6MpU4jNJmNjun79y3eFX1AvSKhWVFQEIDQ0lIgGDhwIQFNTE0CT3ZKp7piamk6ePLlNmzaCX01M4O6Ou3cRGgoNDTg7o2NHeHoiNlbsSqyt4eeH2Fj07o05c9CxI7y9kVl5ZHJiIvr2xfjx8PKqUf2xsbHt2vUdNuzv9u0HHj16QzD5aWqKXbtgbY2YGDx6VLq/ozQKVIoUKXVBWRkBAXj7FpaW6NsXPXsiKanWlVRghBoaGkZGRrt371c/SeMAACAASURBVM7Kytq9e3fXrl319PQACHai0NHRqbNynDx5curUqQsXLnz69GmFBZ48ebJw4cKpU6cGBwfXvbl6jsD/XryAnx8yM9GrF2xs4OuL1NRPZbKzs3NyckTbroYG3N3x/Dl8fXH3Ltq0wezZuHfvU4EPHz5wudxXrzBgAObNg4dHNRWWlCApCTduYM4cvzdvfktL2/nu3SEXl81PnmDxYqSkICwM7u7iXXovRYqUJoieHs6cwe3byMqCqelnYaXFxcXVn1/hOHHv3r1sNhsAm80+fvy44OCmTZt0dHTqPig+cuSInp5eYGCgp6enlpZWSkrKFwWSk5M1NTV/+eWXwMBAPT29o5XN/ja0qdEawuVSWBi5uFDz5uToSP7+NH++h45OD21tu59+Wiu+dpOTS1OGWluTnx8pK7cBWgItZGW3//33p2JZWfToEZ0/T/v20dq1NG8ejRxJ1tZkYEDy8mRsTD17Ups2q1is4wABUY6Ok4QQ02TnnRhEOjVaK5psF61vU6Nfc/o0GRuTrCx16rSTzTbatOlStadUmnT7wYMH/v7+0dHRZUcOHTp06NChuqu0sbHZs2eP4PWoUaPWrVv3RYG1a9eOHj1a8Prvv/+2tbWtrKpGaYRlZGfT/v00YEAyi2Uv2ARRXb1vdHSqWDshl0shIWRu7gf0AvhAEdDG2TnLwYE6dCAVFVJTIwsLcnQkV1datYq2baN//6XISHr79lMQTXJycqtWtnp6o4yMusbExAgho8leZRhEaoS1osl20fpvhAJ8fQmwBD789FP10X2Vpmu0tLS0/HxTn2+//bYOI9dSSkpKoqKi+vfvL/jV3t7+8uXLX5S5deuWo6Oj4HX//v1nzZrF5XJlm15mSTU1uLigf/9ia2sVQbLZvDwVR8finByYmaFVqy9/KtyP8NSpM+HhkSNG2A8Y0P/rd4uLkZCAV6++/PnwIQowAViAPKBoZfW2Z89mhoYwMYGqavXK9fX1X7y4mZycrKen1wT/cFKkSGGcRYuwZEk2n9/M0LD6xYaVXqRev369e/fuR48eFRQUnDp1CsCpU6fU1NTKPEw4UlNT+Xy+lpaW4FcdHZ3k5OQvyqSkpJQV0NbW5vP57969M6rkyVJkZOS4ceMEr3V1dTdu3FgXefUQTU1NOzuN//6bAPB799YPDGxeVJSbnMx+9YoVH89+9Yp19So7Pp79/DlLTg5mZmRqyjc15Qte3L69f8uWC9nZrn//vW7TprS2bYfFx7MFZwn+TUlh6+vz9fXJwIBMTfl9+tCUKXwzM2refHWLFt2BWUC6jMyHxYuNgVyBntzcmipv3rx5YWGhcJ+aiESbE05KtXA/wrSQhkGT7aL5+fl8Pl+ILTskj5WVeXS0Y0mJL9Cp6pIVG2FkZOTgwYNZLJaJiUlZmGhMTIy/v/+TJ0/qokxJSQnlnl4WFhYqfzWQUVJSKisgCF79ukwZRkZGZUNVVVXVKko2XIKD/WNiYlgsVseOHQEoK0NDAxYWn5UhordvEReHV69k4uJk/vsPhw+zIiPPlpRsBYw/fDCaP9/P2nq8qSmZmaFPH0yZQmZmMDbmy8gAYAGszyOnDNPTH7m4uGhra+/dmyDJD1v2cYioUf416y0CF1QUYbqgRk2T7aKCT90gjPDu3bNBQUE2NtV36YqNcM6cOV27dj1x4sS9e/emfNwe2MnJacWKFampqbq6ukIr09DQUFZWTkxMFESfJiYmfj3UMzIySkxMFLxOTExUUVFRV1evrEJDQ8MJEyYIraehYGVlVW0ZExOYmKD8iH3u3HZ//32mpMRNUfHM6tXtV60CUNZ9q+nHmpqax48fl5WVZTOx1wMRsdlsRppusrA/wrSQhkGT7aKCT90gjBDAhAkTeDxetcUq+CtmZmZGRUWtXbu2efPm5T9ty5YtAdQx4yiLxRo3blxAQACAgoKCo0ePjh8/HkB+fv6hQ4cKCgoAjB8//ujRo4LXAQEB48ePbyj/6fWN9etXjRp1p0WL7t9+m/7TT/OZliNFihQp9ZEKRoSC2Ui1r7YVyMzMBFD32AcPDw8HB4cHDx68ffvW3NxcsL9zenr65MmTX79+bWRkNHbsWH9//27duhkYGMTGxl66dKmOLTZZ1NTUjh3bxbQKKVKkSKnXVOBqurq6Ojo6Z86csbKyKj8UO3z4sIqKirm5eR2bbNOmzbNnz+7cuaOmptalSxfBQUNDw/j4eH19fQCysrIhISHR0dE5OTndu3eXPrSQIkWKFCniowIjZLPZS5Ys8fT05PF4RkZGfD7/8ePHR44c8fLyWrJkiYKCQt1bVVJS6tu3b/kjMjIygqlXASwWq8wjpUiRIkWKFPFR8Twnh8NJS0vz9PQUxFJ37NiRxWJNnTr1119/law8KVKkSJEiRbxUHPLEZrM3bdoUGxu7b98+Ly+vP//889GjR3v37pWXl5ewvmo5efL2jh07mFYhRUrDZseOHQoKLZSV28yePZtpLVKkSBoWCb1Jaz0gMDDQxUUXmHXrVpCdnR3TchobhYWFsrKyjKSGIaK8vDzVmqSxkfKRU6fORETcHTPGsWfPnrU9l8UyAP4A1IGZ58/vGzRokDgUNiaabBfNy8trKOsIAfD5fB6PJ1fdNuIVX+MSEhIqW3vRqlWrukoTMYOByV5e+06ckBqhlKaLr+8uD49LWVmTdu9edfDgyiFDHHNzkZuLnBxkZyMrq/R1bi4+fEB29qdfMzMRH38JMAImAQDGjB8f+/33g+ztMWhQxXn7pEhpZFRshHZ2doLteb+m/o0gE4ELJ08eNTKCmxtWrED9m76VIkWMpKfj3j1s3PhvVtYeQC8jQ2/UqAMlJY7KylBVhZoamjdHs2ZQUyv9VV0damrQ1YWaGhQVERWFuLi+wBvgPqAOXNbXn3f0KHbsQHExFBWhqwtzc9jYwNER/fujwgmCjIyM2NhYa2vrJrjAXEojoGIj/Ouvv8pniczJybly5cq///7722+/SUpYTWGzB40e3VNDw/TIEWzYgLVr0akT1qzBx/yjUqQ0Nt68wb17iIoq/TcrC126oHlz8+Tk8zyei6Li2eXLzX/+GVVbUmoq9u7Ftm1o1QpstpypaZ/4+JEA1NV7Gxp2PHUKqqrIz0dYGK5cQVQUAgOxeTOKiqCgAB0dtG6Nbt3Qvz9GjMDvv3t7eu4ATOXl4xMSIuuSeUqKFEaoxTNCLy+vs2fPhoeHi1NP7QgMDAwNDT1w4ACAqCjMm4f8fLBYiImBoiLGjcP69TAwYFplg0X6jFDCxMbG7tsXZGpq6Oo6ufxTjbdvcfdu6c+dOygsRMeOsLYu/enQAWw2srKypkxxj46OcXTs5+fnXcVDkefPsX07AgLg5IRvvoG7O6ZOhYcHeDwul8uVk1OcNQsvX+L0aXyVVAPv3uHMGURE4MEDxMfjwwfw+SCyBUIBLeCPb7+NPnTIX0z/P/WKptlF0UifEVa6H+HXCPJ/xsbG1vwUcfPFfoQ8Hvn7k74+LVhAv/5KxsbEYpGxMXl7f9oqT0rNKSgoKCkpYaTpJrjZ25s3b/T0urJYR5SVVw8ePCM4mDw8yMmJtLXJwICcnIjDIX9/iokh4TaD4/MpLIycnEhfnzw86P17uniRdHVp//7SAmX7EfL5tGABde9O6enVVxsXRzIybYFsgIAdsrJew4bR+fPCKGxYNMEuKqCh7EcogMfjFRcXV1usFhP6WVlZAOrzziNsNlxd8egRiovh54d16xAbC3t7eHpCXh49e+LatdKS+fn5ISEhlT0HlSJF8oSGhqenuxBNyM9fe/Hig8BAKClh4UI8eYK3bxESAi8vuLqiY0fU9l68sBD796NTJyxbhpEj8eoVPD1x5gwmT8aRI3Bx+bI8i4UtW9C7NwYNwse9ZyrFzAwLF37DZluz2c5ycr95eEx5+xbDhkFRET17YscO8KvfDE6KFKap0B7j4+Njy/H06dPTp0/b2Nioq6sXFRWJ2rOFp4od6q9eJUtLcnCgJ08EJalzZ2KzSVOTxo9/KitrzGaPYrNb7Ny5S6KKGxTSEaEEKCmhU6fo229JTe22vPwooAh40KFDf5FUnpxMHh6ko0NOThQW9um4lxeZmdHjx18o+XKH+tWrqUsXSk2tvqEnT54EBQXl5eWVHTl5kvr3J3l5kpWlzp1p+3ZiqCuJi6bTRb+gUY4IKzZCPT29ry3T0NAwJCRE1DrrRBVGSEQlJeTjQ1paxOGQ4Auenk5z55Ks7BIgCCAgVlW1neTkNjSkRihWYmKIwyF9fbK2Jh8fSk2lDRu2m5ra2tqOePr0qRAVnj17bvXq327dukVEd++SiwtpapKbG5WvrKSE3NzI1pZSUr48/WsjJCIPD+rQgd68EUJOKSdP0rBh1KwZsdnUqhWtXk2N4w/bFLpohTRKI6w4WCYkJKR81KisrKyxsbGlpWV9y39dPlimMpKTweHgxg1s2YLhwwHAzm7o7dvTgYlArJzcxOLiOxKS29CQBsuIgzdvcOwY9u1DRgYmTcL336NtWxFUu3Xr7jVrzn/48K2a2ta2bX/OzHSYPRtubtDQ+FQmJwcTJkBREQcOVLBAsLKNeb29sXcvLl7EVzuH1o6rV+Hjg4sXkZMDMzOMHg0OBzo6/HHjpkZEPLCxMQ8JOVAPc1dVRiPuolXT1INl6iFVjwjLc/kyWViQkxO9ekXR0dGyskZsthOL1UpWNtzYmJ49E7fSBol0RChC8vMpKIicnEhLi1xcKCxMyJiXyrCxGQGkAATc6Nt3IZf7ZYHXr6lzZ1q0qNLAsQpHhAI2bCAzM4qLE43UJ0/IzY309IjFIgWF9cB04Cmw0MFhvGgakAiNr4vWkEY5Imwqq1/t7XH/PhwdYWeH48c7zZv3vZ5e3Nix/ePiuhobo0MHLF/OtEQpjREeDxcuwNUVxsbYvx/OzkhKwv79cHSsdcxL1WhptQYuAFBSujBiRBsZmc/ejY5Gr16YMQO+vtWsL6yQpUuxbBkcHBAbKwKp7dvDzw8pKXj6FGx2GLAYaAf8dPVqyunTIqhfipTa8mnWKzQ01MvLq9oTrly5Ik49YkRODu7uGD0aY8YEPHz4gc+/Gxx8UFl59c2bW/76CwsW4MABXLokmnkqKU2TkJCzc+euLC7mTZky7vvvPQIC4O8PIyO4uGDTJujoiKvdU6dw9+4v3bsvTEnZ3L9/rx9++Oy27tw5TJuGP//E6NHCNzFnDmRk0K8fwsJgYVFXwQLMzeHkZHz0qA+wCtiurNx93DhwuTA1xejRWLIEJiaiaUiKlKr5dHMoJyenWgMY1CoSTE3RvfsjPn8MoMjlOkdFPQQwaxaSkyEYGq5cybREKQ2W2bNXvHlz+f37+1u2PBg69L6iIiIicOcO3N3F6IIHDuD77xEcrH77dkBi4p2AgC3ln4j89RdmzkRwcJ1cUMCsWdiwAYMHIyamrlWVcfjwbienQnX1MQMHJrx/711UhEuX0KsXDh5EixZo3hzDhyMwULoGQ4qYkcAsrfio+TPC8pw/f0FDYxBwkcWa4+q6vvxbf/5JcnJkbEwvX4pOZYNF+oywViQlcZWVuwMEkKoq5+zZcxJodNs2MjKiBw8qeIvPJw8PatuWXryoUVVVPCMsz+HDZGhI9+/XUmjtSU8nb2/q0YPk5UsjTpcsocREsbdbQxpiFxUJ0meEjYRBgwYGBXFmzjy7fHnf8+d/fPbs01uzZyMlBcbGaNcOP//MnEQpDYp377B8OSwtZQwNuzdv/r2Kyi/6+lf69u0j7na9vbFtG27cgKXll28VFmLSJFy4gBs30KaNKBudOBG+vhg8GLdvi7Lar9HUxLJluHkTVQwTu3cfJivbUkmp1cmTJ8WrRkqjpqrI+NTU1Li4uNzc3PIHHR0dxSxJEjg6DnR0HAigbVsMH46bN1GWKFhTEzdvYscOLF6Mgwdx7Zo0W6mUSklMxKZNOHgQkyfj8WMYGGwPDw/PyMgYOvQnZXHuYESEH35ARASuXq1g0jU9HWPGwNAQFy5AHCuevvkGSkoYNQonTqD2Wx8KQ//+6N8f/2fvzuNizv84gL9muu9Ld+lyFlHJFUI5onLmFsvKrWWXLJawCMvKnWsduW/JWfySMywiYXOnS3TfM/P+/TFYm9I1M9+Z+j4f/sj4zvf7Kt969/l8PweA9HSsW4djxzBmDHx9LxBpAy/5/H+GDPEoLOwriSisWqnMdmJqaqqbm1vlj2dK9bpGS5kzh5yd6as1MT5JTKSmTUlenlaurOEVZBXbNfodL1/StGlUrx5Nm0bJyRK9dFERDRlCnTtTVta/LxYXFy9YsKJr1yHz529p2JCmTavy9IxKdo1+cfYs1atHly5V7Soi1L37VGCtsC8acNDVJScn+vFH2rfvP18ZMZH+W1RM6lDXqJ+fX1xc3N69ez09PceMGXPu3Llp06Zpa2vv2bNHnEWZGb//jsaNMXp06QfypqZ4/BiLFuHXX2FrC3ZdUpbQy5cYPx6OjlBRwdOnCA6GkZHkrp6fjz59kJ+PM2egqfnv63PnBq1YkX3p0qLFi6+7uOwLDhbx9Ixv9eyJI0c+dcAyIijoBy53BbCNw/nBykpz8WLUr49Ll/Djj9DWhqIijIzQuTN++w1RUWUPtzlx4kRYWJjEg7OkT5klVEVFZe/evUQ0evToOXPmCF9fs2aNnZ2dVP0uIJIWIREVFVHnzvT5Ey3tS9Nw9WoioqysrNTU1JpfVPqxLcJSEhLIz48MDGjBAvr4kYEAHz+Siwv5+paxbqejoweQBhBwc9CgSdU4eVVbhELR0aSr+9jc3NXIyHHYsEl8ye7zEhMT06/fkHnz5pW6bl4enThB06eTqysZGhKXSxwOaWpS06bUrx9t2EDv35ORkT2X25PL7WZm5lSNS0vnLSoBtbJFWEYhTElJASBc7XD8+PFTpkwRvv7hwwcAz549E23QmhBVISSi9HRq1Ig2bSr3gMWLSU6O1NV/43ItudxGdnZdRHJdacYWwi8ePaKRI8nQkBYsoMxMZjIkJ39aGqbMn0L+/os5nDnAU01N31279lXj/NUrhETUuHEPIA4gNbWfDxw4WI0zSMCDBxQURL17k7U1qagQ8BLoJexW5XBce/V64u9PS5bQ1q0UHk4PHlBGxvfO5uMzHDABzBUVLST0CUiNWlkIyxgso62tzeFwhJsumZqaXrlyRfh6QUEBpHsbpprQ08OZM+jUCTY26NatjAPmzYOnZ4aDw3HgOcCNj+8UGRlZ3pNUVq3x8CFWrsTFixg/Hs+e/ac3UpJevkT37hg+HIGBZR+QkhLQunWwunpg//5uvr5DJZmtsDADaAwgP7/ly5fvJHnpyrO3h709AgI+/fXmzcz27dOECy0T5T19qvHwIfLzUVSE4mLweJ+6UrlcyMtDQQHKylBVhbo6tLSgqYkLFy4DdwGj4uJhnp6ep0+fZuwTY4lCGYVQSUnJzs7u9u3brVu39vDwWLhw4YoVKxwdHVetWqWlpdWw9q68YmODAwcwYAAiI8sYjw5AXz+fy+UIBFwAAoFJfHw2Wwdrmfz8/FGjpt+6dbddu1Y///znH3+oXLmC6dOxeXMZq1RLzKNH6NULv/2GcePKPmDDBsTHK9y69Qsjq+KPGjVw3bphWVkdFRS2Dx58jIEEVde2bcsWLXRjY20BQatWDW/dMvn2mMxMvHmDxEQkJSEtDWlpSE/Hx4/IyAAgBwi36HGIjFQaOxZt2qBNG9jZgYk16lk1VmY7MTQ0dM2aNcKPp02bJlxoXE1N7cCBA6JstdaYCLtGvzhwgKysyh0HaGnZhsvtzuUOlJPryOXy+/Wrbbusfa0Odo3+/HOgomIwQFxusKbmwg0bqLBQ8in+IyqKDAzo0KFyD4iNJX19qtbGTf9R7a5RIoqKitq2bUezZon799c0hiS9efMmMTGxGm9UUbECvIFFgKmDw9PNm2naNHJyImVlcnKiadNo1y569EjkeaVCrewardR0iPfv38fExGQy9WykfOIohEQUGEhOTpSbW/a/Hj16dPfu3Xw+/9Ah0tAgTU06flzkEaRCHSyErVuPAOIBAuLd3IZLPkApYWGkr0/ny1+jJieHmjSh0FARXKsmhVAoJoYMDSu1kW8tMGjQoBYtWkRH35o2jRo1+rS4T04ORUfTmjU0ciRZWZGREXl60oIFdOoUpaf/+94nT55MmBAQGLgiSwLzPEStDhXC6v2WJHliKoQCAY0aRZ6e9O1eNqXw+eTnR1wutWhBSUkiD8KwOlUI798nHx/S0zusotIbOKGj03v//iOSDPCt0FAyNKTr1793zMiR5OcnmsvVvBASkb8/jR4tkjjS7utbNDSU6tWjHTtKH/PuHZ06RQsWkKcnaWuTtTWNHElLlrzX1XUEzigobG7dupekc9dYHSqEhoaGTk5OISEh2dnZog4mSmIqhERUXExubvTzz5U6OCGBmjYlOTmR/UiSEnWkED54QD4+VL8+rVlDBQUUGXlp1qxFkZHMTRQnIqJ168pdRPSLrVupWbMy1oKoHpEUwrw8sramCxdEkkiqlbpF4+PJ1pb8/KioqOzjS0ro/n3asoV69DjD5f4mHLBar167yvyYlip1qBCGhIS0atVK+FzQ19f30qVLEp4eVEniK4RElJlJzZrR+vWVPX7tWlJSIkNDio4WUyJJq/WFMDaWfHzI0JCCgqjGJUAEkpOTN2zYfOrUqWXLBE2a0OvX3zv40SPS16fHj0V2dZEUQiI6c4asrMp9slBrfHuLZmeTjw+1akUvX37vjS9evKhXrwOQBtyVk2u3YUPFPU9SpQ4VQqH4+PgFCxZYWFgAMDMzCwgI+KeS69hLilgLIRG9eEHGxnTqVGWPz8ujfv2IyyVXV5KmWXDVVIsL4cOH/5bA/HzxXacKUlNTTUxaysltUlScrKc37ftP2goKyN6edu4UZQBRFUIiGjqUZs0SyZmkV5m3qEBAa9aQsfH3HusS0ZEjJ5s169q5s8/Zs0969KAWLejKFTFGFa06VwiFeDzemTNnhgwZoqyszOFwRJFNZMRdCOnz8/9796rwlqgo0tcnRUWZX6S0VhbCL1PjpacECu3fv19J6Q9hj5mhYQVrnYwZQ0OHijiACAvh+/dkZER374rkZFLqO7doVBSZmFBAAFWyH+3UKbK0JE/PCvoApEStLIQVb8MkJydnbm5uamqqqalJwgmodYmzM9avh5cX3r6t7Fs6dUJaGmbNwq+/wsYGT56IMx+r0h4/hq8v3NxgZ4cXLxAQABUVpjN9xcjIXCC4DQiAl5qa35sPePAgoqMREiKxaFVWrx6WLcOYMSgpYToKEzp1wr17uH0b3t7CSYcV8PJCXBycnODoiMBAFBWJPyLrv75XCNPT09euXevk5NS8efOtW7d6eXlFR0dLLJn0GDgQkyejTx/8d0OqCixejORkmJvDzg79+4PHA4/Hu3//fnFxsdiSssoWHw9fX3TtCmtrPHuGgAAmZ8eXqagIa9a4WFg0NzZ2btr0hwMH1pV3ZEICpk7FoUPQ0JBkwCobPRoGBli7lukcDDEwwPnzaNYMDg6V2rhRVRWBgbh1C48fo3lznDkj/oisr5XZTgwLC+vbt6+ioiKXy3Vzc9uzZ0+eqIamiZQEuka/mDiRPDyqM33+8GHS0CBl5edcriWX20FOzvSCjAyqk9Gu0Y8fP+7Y8deJEyf4fP6LF/8uky21U7by8qhbNxo8mCrswiksJAcHCgkRSwwRdo0KvXxJ+vqUkCDCU0qRSt6iJ06QoWHV/ssiIsjWljw96fnz6scTn1rZNVru9In69esHBAQ8F9t/RXJy8tOnT2v4BZVkISwupm7dyM+PkpKSXr16VaX3lpSQhkY/4AxAwDUjI0cxhRQtWSyE2dnZFhZOCgp/qqlNt7AYLSyB0rcUxL9yc6lrVxo2rFK/Y02cSD4+4koi8kJIRMuXU9euVd4ZUSZU/hZ9+pSaN6eRI6sw0aW4mNasIT09mjZN6obd1cpCWHbXaHh4+KtXr4KCgqytrUXeBhUIBD/++GOLFi369Onj4OAg3OyiVCPVw8PDwMCAw+Fcu3ZN5AGqR0EBR47g+PGFjRoNc3aePGDAj5V/r7w8dHSSAXUAgFp6uuHNm2KKWdddv349Pb1HSclPeXmr09Pj4+N5gYHQ0mI6VjmystCtG2xssGdPxWtUHjmCixexdatEkonIjBnIykJoKNM5GNWoEW7dgrw8XFzw/Hml3qKgAH9/xMYiIwNNm2L3bjFHrPPKLoROTk4csW3ree7cucjIyKdPn8bHx9vb2y9atKjUARwOx9fX99q1azo6OmLKUD1cbi6ffzY399L796ejo7MfP35c+feuXz+Hyx3G4Yzhcvs1aDDVxQU2Njh5Unxh65z8fGzbhunTTQsLHwB8IFVbu0RXV3qXQM7MRI8eaNkSISHgVjRq7cULTJ6Mgwelt6iXSV4e27dj5kykpTEdhVEqKtixAxMnom1bHKv0suQmJti9G3v3YtUqdO2KR4/EGbFuK/f7Lzw8fNy4cb169er2XzW/5IEDB4YMGaKtrQ1gwoQJBw4c+PaYoUOHNmzYUHzFuHo4HI6c3KeP+fyqDaD18vJKTIwJCWn36lVUfLxHYiJatMCAATAywubNoo9ap/zzD2bPhqUlDh7E8uXN5s93MzZ2trbuFxq6hulo5Xr/Hp07o1MnbNxY8VbyJSUYPhzz5sHRUSLhRKpFC4wYgZ9/ZjqHFPDzw7lz+OUXODnN1ddvYWLS8tixivdvEo5BHT0a3brB3x/Z2UhISNi3b19CQoIEMtcRZRfCGTNmeHp6nj17Vhy7D7569crGxkb4sY2NTUZGRnZ2drXPlpmZefezR2L+lUlNTW3kyF76+l01NHqXlOiYmdlW6e3Gxsbjxo0zNzcHYGyMY8eQmYk+feDvDy0t/Pbbpy3QWJUkIPsONQAAIABJREFUECAiAl5e6NQJAGJicPEivLwwf/70pKS/nz+/3rlzR6Yzli01FW5u6NkTK1ZU6vhZs1CvHqZMEXMssVm0CDdugN2zD4CTE7Zu/fvhw2fp6Q+Sk6MnT/6tMu/icuHri7g4ALCxuezgMHLMmJR27XwjIi6JN26dUUbHEZ/P37x586RJk9auXSv3pQVUFY8ePVpb1rjpZcuW6enp5ebmKn/eNk1FRQVATk6OZnU3PH3w4MG4z7u0mZiYlNm+FKGFC2dMmjS8uLg4KKjBmDG8nTsLa3jCVauwfDmCghT//FPhjz84Q4eWBAUVScng/sLCQnl5eXkmNlgjou/8EpaWxtm7V2HbNgU9Pfrhh5IdO0qEMwKrNL+FKSkpHC8vlQEDeLNnF1cm8Pnz8seOKUVH5+fliXcWL+8zcZx87Vq5iROVnJ0L1NRqyVzk79+i35Gf/05V1TIrC4AGjyefmZlZyW8xRUUsWYKYmL9u3gwB7IuKegQFLW3btnU1MtREfn6+QCCQtu668ggEAiUlpQoPK+M/ID09vaCgYMyYMdWrggC0tbXbtm377evCQIaGhhmfZ5l+/PiRy+UaGBhU70IAXF1d9+7dW+23V4O6ujqAkBB07Iht29R/+kkE5wwKQlAQNm7EokUKe/YoeHlhyxbUqyeCM9eE/GeSvzQRcTgc4Zf6a3fvIjgY4eEYOBCnTqFFCw6gBFR8o0uJN2/g4YGJEzFjhiKgWOHxb99iyhQcPYr69dXEnU1YBZXFs7dvr15wdcWKFWqrVonj9Awo7xatkLu7u5FRENGvxcXv5OVdVFS0K/GD+l/29ka3bz/h8+2BJ1yucTUC1BCHw1FVVZWhQsjn8ys+7tuBpDwez9DQ8OTJk6IdxvrFnDlzhn5eHmr//v329vblHamrq3v16tXvnEqS0ye+9fo1GRtTVJSIT7tvH1lbf1qwlNk5WNIzfSI7m0JCyN6eGjWioCD6+JGRUDX14gVZW1NwcGWPLykhFxfJLdQnjukTX0tPJyOjCnaVkiE1mepaUFBw+vTp6OhrAwfSkCGVXYlNKD093dnZw8DAydbWw8Li/bRpFU8/Fa1aOX2i7HmE+/bts7e3f/n9RdSr6+XLl5qamlu3bo2MjLS2tt62bZvwdS8vr9DPG4weOnQoJCRETU1t1qxZISEh5W1fyWwhJKKICDIxIXHs3njpErVo8WmnQxeXsVyuGZdr2qePRLd6Y6oQnjhxSkXFRk7OrEePwU+eUEAA6emRpyddvCjDM9KePCFzc9q+vQpvCQggDw/JfcriLoREFBpK9vaS/sEtJiJZDreoiNzdadKkar49PZ169aIOHcTyI6g8tbIQlj1Y5vjx4ykpKY0bN27ZsqXIR41aWlqeOXPm3Llzy5cvDwgIGDt2rPB1Z2dnU1NT4cdxcXF3794dPny4cCxMYWFNH8WJiZsbJkyAjw9Evm5aly64fx8xMeDxkq5dixMIXgkEr8LCbrx+/VrEV5I+gwZNKygI4/NfnD+f367dNVVVPHyIsDC4u1c8wFI6xcfD3R2LF2PMmMq+5dw57NuH3btl9VMu0/DhqF8fq1cznUNqKCri6FHcuIFly6rzdj09nD6NgQPh7Izz50Udri4p9/GPvb29+K7q4uLi4uJS6sXffvt3AFVgYKD4ri5a8+bh3j3MmoU1Yhir7+SEDRv+6dpVQyCQAyAQNGnfPr9/f0ybhoYNRX85ZvF4iIlBZKRwpeamAIBOQ4bsnz+/9K0iW+7dg4cH/vgDI0ZU9i2pqRg7Fnv2MP+cWOTWr4eTE/r0QZMmTEeRDpqaOHcOHTvCwACfGwVVwOHA3x8tW2LYMIwciSVLUN2hHXWbBBqn4sN416hQRgY1aEC7donl5Hw+X0PDhsMZx+GMVVdv6udHZmbE4ZCqKrVtSxs2VGf500qSTNfo8+cUEkI+PqStTdbW5OdHurqOwAxgD5drfq9KO2BJnzt3yMiIjh2rwlv4fHJzo99/F1umckiga1Ro9WpydZXhXm4h0e4UlpBAJiZ09Gj1z5CaSu7u1KULJSeLKlTZ6lDXKKtKtLVx7BhmzhTL0g9cLjc9/fGyZTYrVjTJyIgNCcHbt8jPx8qVADBzJpSUYG6OkSMRGyv6q4tJSgoOH8b48ahfH9264epVuLvj0SM8f46QELx7d83Pr6BLl0PR0QdatmzJdNjqu3oVHh4ICUG/fpU6ft685WZmrYyMumZnx86eLeZwzPH3R3Exdu5kOoc0Ea4zNWECqr3Bj4EBzp1Dp05o3RpSszCl7CivQt64cWPw4MHNmjWztbUVvrJmzZrtVXrWL35S0iIU2ruXGjZkYH3nBw/oO83E7du39+s36Gi1ftWsSYswLi5uw4ZNt27d+vJKdjZdvEgBAeTkRPr65ONDISH04kXZbxf3DvUSEBVFhoZU+Y1Grl27pqPjA/CAf6ytO4gzWtkk1iIkogcPSF+f3r2TzNXEQhy36KVLZGBQtW3Av3X69KftVsTUbKuVLcKyC+GpU6fk5eWbN2/u7e1tamoqfHHTpk3m5uZS9SWQqkJIRJMnk5cXY30+BQW0YQO1bUuqqsTlkpkZ2dou43BcgL84nOZ//LGm6iesZiG8cuWqnl5bDmebjk6PWbP2C4ufhga5u1NQEN25U/F4cRkthDwe79ixY1u2bDt8+GO9enTpUhXeu2/fvi871BsZVbBDvThIshAS0ezZNHiwxK4memK6RffvJzMzquGA/TdvqF078vYWy0SjOlQIra2tR4wYwefzL1++/KUQPn36FECiJAfqVkTaCmFxMXXsSEuXMp2DKDqahgwhDqc78Aog4LaCQv9evWj8eAoKohMn6NmzCs6wZ88eDQ1LPT2b27dvf//IkhJ6947+/pvCw2nnTlq2jJo0mQZcBQhI0tDwDAyk6OiqPcuU0ULYu/dIdfWf5eSC5eScIiKyq/Te169TFBRacjjb1NX9hw6t7oD6GpBwISwspCZNSGzTlcVOfLfoxo3UoAGlpNToJCUlFBBAFhb0VaeMaMhQIczIyOjbd/ilSxX8BCOiMkaNpqWlvXjx4vDhw1wu9+vlA4RzG1JSUr5McmCVoqCAQ4fg7AxHR/TowWSSDh3QoQNu3855/vwmYAFc09TU+/ABT54gMxN5eSguBhHk5aGiAnV16OujXj1YWKB+fTRtCiWl+JEjZwKrgAxnZ6/k5OS0NCQnIzUVaWlISsL79//+NSMD+vowMICJCfT1YWwMCwuz58/vlJS4cLkxvXqZL1jA5JdCYvh8/u3bT3NzdwNQU0vl8a4DVbgJFiww9PYOd3E5ZmHh2rdvH7HFlBZKSti8GaNGoUsXaGgwnUaaTJyIxER4eeHSJVR73Rh5eQQFoW1beHvj11/h7y/SiDKiQQOXDx+GubmZdOlSwZFlFEJh8SMqvSRgUlISAFUpWQdTWhkZ4eBB9O+PGzdgZcVwmLNndzk7e+fmBujqaj5+fKnUWPzsbDx6hIcP8eoVnj9HSgquXkVmJnJzUVT0F+ABDAMAnLSzSzMzMzA2hqEhDAxgbg5HRxgZwcgIBgbQ1y89162gYMqAAePu329lZWW6bt02iX2+zCKSy8/nAamAnrLyfVPToZV/78qVePAAV6+aqKrK7NLaVefqCjc3zJyZPWpUXJMmTaRt2zUG/f470tLQty/Cw1GlBdhK6dsXdnYYOBC3byMkBGpiX6dPily5go8fCZj77l0ldjMos51oZWXl7+9PRFFRUV+6Rv39/fX19Xk8nuharjUlbV2jX6xcSQ4OlJ/PdI7qOnz4MNAI+AC8BEyDgnhFRRINIHNdo+np5O5OrVtfsbRsZ2Li+PvvVXgie+ECmZrSmzfiS1cpEu4aFbp585GcXEstrSmGho43b4q6F0+cxH2L8ng0YAANHVq1BdjKVFBAfn7UpAk9fCiCYNLcNZqaSvPnk709KSoSl0scjiNwb+3aimtW2YVw9+7dAIYPHx4YGGhgYHDkyBEfHx8AwZVfJ1EipLYQCgQ0aBCNH890jhpo3borYAKYDhrk7+NDDRrQoUOSu7psFcIHD8jamgICqBq/Jb54QUZGol+xthoYKYRDhkwGogACYrt1Gy7hq9eEBG7R/Hzq2JEmTxbN2XbtIj09+vPP1AkTAvr2/TEmJqZ655HCQvi//1G/fqSrSwBpapKHB+3ZQ3w+nT17VkOj0ebNZyo8Q7nTJ7Zu3aqvr/+l4aiurh4UFCRtn7/UFkIiyskhW1v6vJCqTPp61GhkJNnbU9eu9OCBJC4tQ4Xw4EHS06Pdu6vz3pwcataMQkJEnalaGCmEY8f+zOGcBgi45u09RsJXrwnJ3KKZmdSyJQUFieZsjx+TsrI7h3MIuGVg4PyuWvNXJF8I+Xx+TExMamrq1y8+fUrTp1PTpiQnRwoK1LQpTZ9expqrNRo1KlRUVHTjxo1jx47973//k86fStJcCIno6VMyMKCKBl1Kr1LTJ/h82rWLjIxo5Ej67z0pejJRCAUCWrCA6tenO3eq+faBA2niRFHHqi5GCuHbt2+trJz19XvLyTkdO/aPhK9eExK7Rd+9Iysr0fxKzePx9PXbCKfoqKoGnjp1qhonkXAh/PDhg4qKBZfbicu1mD178Z495OFBmprE4ZCuLvXrR9HR33u7CAqh9JPyQkhEJ06QhQW9f890jmopcx7hx48UEECGhhQURIWF4rq09BfCrCzq04c6dar+7wS//UYuLiThh6/fwUghJCKBQJCcnLxrF9/eXoq+GhWS5C36zz9kYlK1hfrKY2vryuGcB55wue0mTnxdjW9hCRfCsWP9gN8BAjKBNurq1LEjbdpElbxVK1kIy150+8qVK8Vl7aegqalpaWlZk31065o+fXD1KoYOxblztWQxXB0dBAVhzBj8/DO2b8eSJfDxYTqTxD17hr590bEjDh+GgkJ1znDiBHbvRkwMFCvenbeW43A4RkZGvr44cQJLl0J21tuXnAYNcPIkevWCvj46dKjRqc6d2/PLL0tTUz9Mnbp03776Tk746y84O4soqKjdvImTJ80BFQCAvJxcak6OeK5UZnk0NDT8zltcXFwSmN0x9jPpbxESEY9H3brR/PlM56i6CleWuXiRmjUjNzeKjRXxpaW5RXj6NBka0o4d1T/D48dkaEjVHawgLky1CL9ISiIDg2r2M0ue5G/RyEgyMKD790V5zkOHyMiIAgKq0LsjmRbhnj1kbU0cDtnZpSsoWHK5fbncRhMn/lzV89Soa/TYsWN6enqTJk06e/bs7du3T506NWzYMFNT05MnT27cuNHExKRp06bSMI9CJgohEaWmkrk57dz5/vLly6nifrwmOpVZYq2khEJCyNiY/PwoLU1kl5bOQigQUFAQmZvXaLWODx+oQQPau1d0sUSE8UJIRDt3kqx0kDJyi+7fT+bm9PffaREREaL6SZKaSgMHkq0t3bxZqePFWgjz8mj6dNLQIAUF8vCgV6+IiAoKCo4ePfr06dNqnLBGhbBdu3YrVqwo9eLEiROFVef27dsAor//jFIiZKUQEtHu3fe4XActrZ8NDJyuXr3GdJxKqfxaox8+0LRppK9PQUGi+SkmhYUwJ4f69ycXlxptc8PjkYcHzZoluliiIw2FkIj69qUFC5gOUQlM3aIzZ8bIyTlqac0yMBDlzMtDh8jQkPz8KDe3giPFVAjfvKF+/UhOjjQ1afp0kf0yVP1tmD5+/Hjjxg1vb+9Sr3t7e4eHhwNo1aqVsbHxy5cvRddBW/uFhW0RCNZnZf2RlrYrMHAD03FETFcXwcGIjsaVK7C3x+nTtGrVxg4dBsydu7TMh80yJyEBbduiXj1cugQjo+qfZ+ZM8HhYulR0yWqdjRuxaRPu3mU6h7R6+jSEzw/JylqelrZ98eJNojqtjw/i4wGgRQtERYnqrJVy4QJatoSFBR48wKFDyMrC6tWSfnZeRiEkIgAJCQmlXhc+FxR+rKioqKysLO5wtYm2thqHkw4ASNfUrJ0rHTVujPBwrF+P8eNDZ8++d+3an3/+WTJ3bhDTuWrq/Hm0bw8/P4SE1Oj7MzQUJ09i//5aMmxKTIyNsXw5xo5FSQnTUaSSlpYah/MeAJD28qW6CAeP6OggJARr1mDECIwfj7w8kZ25TDweVqyAkRE8PKCtjSdP8Pw5+vcX70XLVWY7sU2bNhYWFlGfl7sQCAQnT57U1tYeMmQIEaWlpcnJyVW4KYEEyFDXaFJSUsOG7fT1uygoOK9c+ZLpOJVS7W2YBg2aDNwECEhzcPCoxhmkp2t0zRoyNKTLl2t6nr//JgMDevRIBJHEREq6RoX69KGFC5kO8V1M3aKJiYkNGrQ1MHC1tGw7fPhbU1MKCanOkkbfkZFBfn5kbV32PmLV6xpNSkpq1aq7gYHj0qUrUlJoxAhSUiJVVRoxgjIyRJC5PDV6Rvj06VNra2sAGhoa1tbWKioqABwdHYWPZ69duzZx4kR2sEw1ZGVlPXkigr03JaPahXD37v2amiOAODm5Wfr6S6sxtEQaCmFBAfn6koMDvX5d01OlpJC5OR0+LIpYYiNVhfDdO9LXp7t3mc5RPmZv0S+XvnOHOnQgR8cK5pVXQ3g4mZuTnx+V+iyrVwjr1WsGbATuAW05nL8NDWn5cpFF/Y4azSNs1KjRw4cPjxw5Ehsbm5KSYm5u7uTk1LdvX3l5eQDt27dv3769JJuttYampqamJoKDMXgwbt+GpibTgcRj5MghhYVFR48Gubo6WVlN9vZGnz5YvVoGFr8vLi72959/+fL1Nm1c4uIWNmmieO0aVFRqdM6SEgwahHHjMHCgiFLWASYmn6ar3r5dzZmatZv65/2ZnJwQHY2wMIwYgebNsW4dLC1Fc4levRAbi4AA2Ntj2zZ07Vqdk/zzD44fR3Q00tMBTAQA/Ojmtv7ixe2iSSkqEqjJ4iNzLcIv/Pxo0CCmQ1Sk2i3CUj5+JD8/srGpwo7tTP26PW9ekLLyYiCPw1nUo4dofmX186O+fUnKluktg1S1CIU8PGjRIqZDlEMaOi2+lpdHQUFkYEABAaXbcDV05synpmF2NlFFLcK//6Z588jVlczMSEmJOBzS1KSmTUlBoQOwH3jB4bTds2ePKPN9V/VHjbIkIDgYCQnYuJHpHBIhfAi/di1Gj8b48RDX2hCiEBn5sLDQB1AlGsznP6j5CTduRHQ0du0qvWUjqzK2bcOGDXj4kOkcskBVFQEBuHsXSUlo2hS7d+ObLWWrycPj039BixZQVLRRV2/C5ZoOGjQCQGEhTp7E+PFo1w5GRpCTQ6tW2LgRmZno1QsHD6K4GFlZePwYcXE7GjZcr6npOWNGlxEjRogmmQh9KYkHDhwwNjZet24dEdnb2xuXQ4y1u+pkt0VIRAkJ0r4kt6hahF9kZpKfH5ma0smTFRwp+V+3Hz6kkSNJU3OXktJQIEpLa8jWrdXaUeIrV6+SkRH9IyNLSUthi5CItm2jli2pEr/TS5q0tQi/duUKOThQx44ifsjavv0wYCBAQDLQXEWFAFJSIktL8vSklSvp2TNRXk4kqvyM0MrKysfHp0mTJgC8vb2zs7OZq851go0Ntm3DwIG4exd6ekynkQgtLYSE4MoV/PgjQkOxaZNUfOKPHmHFCly4gAkT8Pat7/nz6uHhJzw9Bw0c2K8mp33zBj4+2LEDDRqIKmldNHYsjh7FihWYO5fpKLKjY0fcuYMdO+Dpid698fvv+O6imRV49Qo3buDGDdy/Hw8MAAAYAUUhIfDygra2iEIzikOiaj8zITQ09OzZs3v37mU6SPXNmIF//sGpU9LYdVZYWCgvLy8cISVa+flYtAh79mDtWgwYUMYBRJSXl/dlRICYPHyIxYtx5QqmT8fUqVBVFdmZCwrQqROGD8dPP4nsnOLG4/F4PJ4Uzg9+9w5OTrh4Ec2bMx3lK5K5RWsoLw8rV2L9ekyZgl9/hZJSpd5VUoLYWFy9irt3ER2NwkK0agUnJ2hpXZkxYzAwDbijpHSvsPCFmOOLgEAg4PP5ChUNuGKfETJs+XJkZGDVKqZzSJaqKoKCcPgw5s3DoEFIS5N0gNhYDBqEbt3g5ISXLxEQIJoqSERnzpzZu3evr29OkyayVAWlmakpFi/GqFHsFPsqU1NDYCBu3sTjx2jWDIcPg4guXbp06NChvP9OmE9KQlgYZs9Ghw7Q0YGvLx4/hrs7Tp9GcjLCwhAYiOnTO127dtTObu+AARyZqIKVV26L8OTJk3/88UdcXJyqqmpiYiKAFStWKCoq/iRN39y1oEUIIDERrVvj4EF07Mh0lP8SX4vwq0sgKAibNmHxYvj5/fu6+H7dfvAAS5Z8agVOm1bTqRGlDBjwY0SEcn6+GZd7LDHxsr6+1M8X+YrUtggBEKFXL3TsiDlzmI7ymUy0CL8WGYnp05GaOrGggM/j1Tc2Dtu169Ldu2qlmn0dOsDF5XvfF3l5eaqqqhwp7MIqSyVbhGVPn9i5cycAd3f30aNHm5qaCl/ctWuXgYGBNMyj/0KmB8t8LSKCTE1rtJqzOIh8sEx57t8nR0fq1YvevPn0ijhGIty4QZ6eVL8+rVlT2V09q4TH4xkYtBZu/62h8ev58+dFfw1xks7BMl8kJpKhIT18yHSOz6R5sEx5Cgt56urOwlsUmNOo0dmpU2nfPnr5sgonkfDGvDVU/ekTRDRnzhx/f/+LFy+OHj36y+suLi5paWnv3r0TYblmCbm5YfRoDB8OPp/pKExo0QK3bqFTJ7RqhS1bQIRHjx6Fh4fzeDyRnP/GDXh5YfBguLvj6VP4+0MczR45OTmBoAT4AAgUFR8Z1WRxbtY3TE2xaBF8fdkO0upTUpLT1CwBMgDS1X1y4IDh2rUYOlRkc/BlVxmFMDU1NSkp6Ycffij1uvAbO03yz3PqhkWLICeH339nOgdD5OUREIBLl7B9OwwMfmzRYsiwYXu0tBpWY/Ty8+fPf/klcPny4Ly8vOvX4eWFoUPFWwKFbt4Ej/eHoaGHkZHThAku9vb24rpSXTVuHAwM6twDddHatm2FuXl3AwOHsWMdHBwcmI4jLcp4/KOoqAigoKCg1OuvXr0CoKWlJf5UdRGXi9BQODmhbVv06MF0GobY2eH6dSgoXCBKIFLMz5/Tp8/OIUOmcTjQ1IScHDQ0IC8PNTUoKkJFBcrKUFaGigoUFT+t35aent6x46CUlLny8onLlg3T0zs5dy6OHRP7Ml2xsejbF3v3du3VK0a8V6rDOByEhKBVK3h6olkzptPIJg+Pbm/edGM6hdQpoxDq6ura2tpu3LixTZs2X56IEtHy5cvNzMwasLOixMbAAKGhGDoUt2/D1JTpNAyRkwOXC/6nPuLi4mLle/dAhKwsCATIyQGPh9xclJQgPx9FRSgsREEBioqQnw95eSgrx+Tl9SbqX1ICov1xcTxlZTEO9hF69gweHli3Dr16iftSdZ2FxacRpDdvsmuQskSm7J8RQUFBffv2fffuna2tbUFBwbp16w4fPhwdHb1nzx5ZGSwko1xdMXUqBg7ElSt19/t81Kj+O3c2B0zU1ZMvXnxQ+YkNJSV4+LBx9+4rPnzIBBLr1YMEquCbN+jeHYsXw8dH3JdiAcD48ThxAqtXIyCA6SisWqO8UTRhYWHCVWaEzMzMJLlSaiXVmlGjXxMIqE8fmjmT6RwSHDX6rbdv30ZGRlbvvaGhhxo37tiunXdcXJxoU30rMZGsrWnTJnFfRxKkfNTo1169onr1GN7cURZHjYpErRw1WsHKMomJiampqZqamg0aNJDCtmDtmEf4rYwMODlh1Sr0q9EiXzUlgXmE5SFZmKT1/j06d8aoUZg1i+kooiDN8wi/tWkTdu3CtWuQk2MmgEzcouJQK+cRVvAzzszMzMzMTHSpPklISFi1atWHDx969OgxZsyYUl/TnJycQ4cOXbt2raioqHXr1n5+fiqinfYs9XR0cOAAvL3RogWsrZlOwypLVhZ69sSgQbWkCsqcCRM+dZDOnMl0FJbsY2CJtZycnI4dO+rq6g4fPnzlypXBwcGlDrh27drJkyfbt2/ft2/fvXv3Dhs2TPIhGde6NebMQf/++Gb0Lot5+fnw9ESHDliwgOkodZVwBOnChZstLNp37jxQOKadxaoeBnq99u7da2Njs2TJEgDKysrjx4+fNm0al/tvSe7Ro0fPnj2FH9vZ2TVv3jw/P19VhCsiy4hp03DrFmbOxPr1TEdhfaWgAJ6eaNgQa9YwHaVuS06+QXT+zZvIN2/uDxgw4e7dc0wnYskqBlqEt2/f7vh5Vc2OHTu+fv06JSXl6wO+7il9+/attrZ2Xesa/SIkBJGRCA1lOgfrs5ISDBoEExNs2yaNG4bUKfHxT4qLuwEqQLvk5A9Mx2HJMLG0CEtKSpKTk7993cjISFFRMSUlpWnTpsJXVFVVVVVVU1JSTExMvj0+Ozt72rRpixYt+s6D2evXr3ft2lX4sbGxcUhIiCg+Aymycye3V6+X27at1NTkBAZObdSokcQuzexgmfz8fMlf9/v4fPz4ozIR1q8vlL50NcX7jOkgldWmTWsdnVHv35tyOPfMzJrn5uZK8urSeYtKQH5+vkAgkKHBMkqV2H1KLD/jnj596u3t/e3rR48edXBwUFFRKSwsFL4iEAiKiorU1MpYpD8/P9/T07Nr166TJk36zrUaNWo08/PjchUVldo3iKtlyyKBYGRU1AZA8Pff416+vFnxSuoiIv+ZZC73NSLicDhS9b9JBD8/ZGYiLAzKylIUTFRka9QoADs7u8uXd23ffkhJyXz79jmvXyva2Unu6lJ4i0oGh8ORuVGjFR4mlp9xzZo1e/Gi3N2qzM3NX79+Lfw4MTGRiL5tDhYUFHh7e9vY2GzYsOH7X/F69eq5u7vXPLPUevXqlaKiLdAOAI/X+NWrVw0bNmQ6VF30yy+Ii8OFC2JcrZQuN1IiAAAgAElEQVRVVXZ2dqtXLwRga4v+/RETA3YJSFY1MPCM0MfHJyws7MOHDwB27tzZs2dPDQ0NABcvXrx37x6A4uJiHx8fXV3drVu3fj2Ipm6ytLRUUIgHooGooqKnluxC8UyYMweXLiE8HHWvASAbRo5Et24YORICAdNRWDKIgTLTvn37Pn36tGzZ0tXVdcuWLUFBQcLXV69effz4cQDHjx8PDw+PiYlp3LixjY2NjY3NlxZkHaSkpBQRsXfAgNDevQ8oKBy4fr2uLrzGnCVLEBaGiAjo6DAdhVW+P/9EdjYWL2Y6B0sGMfD4B0BISMiLFy9SUlIcHR2/PJPYv3+/8OlXv379Pn78+PXxdXzLiyZNmhw5EgIgMhLDhiEmpu4uyS1569dj505cuQI9PaajsL5LQQEHD8LZGc2aYcAAptOwZAozhRCAtbW19X8XTdHW1hZ+oKioKNwKilWKm9unJbn/9z9UYiQUq6Z27cLKlYiKgrEx01FYlWBoiBMn4OGBpk1ha8t0GpbsqOtP4GROQADq18dPPzGdow44dgy//orz59n9u2WJoyNWrkT//sjKYjoKS3awhVDGcDjYvh3R0di+nekotZFAIFi5cn23bsPHjVs/aZIgPBxfbcHCkg2+vujSBaNG4bsbCrBY/2ILoexRV8exY5g7F7dvMx2l1lm9etOiRU8jIuZv3/508ODNDg5MB2JVy9q1yMjAkiVM52DJCLYQyqRGjbB1K3x88P4901Fql/Pnr+fmTgUaE02Ni7vGdBxWNSko4NAhbNmC06eZjsKSBWwhlFVeXhg5EoMHQ3aWxJIBKirtOJwNwHN19Q09e7ZnOg6r+gwNcegQfvgB8fFMR2FJPbYQyrCFC6Gqitmzmc5RKxQVYfRoJCZOCgiwdnWdM2+ezYwZE5kOxaqRtm3xxx/o3x/Z2UxHYUk3xqZPsGqOy8WePWjdGs7OGDyY6TSyLD0dAwbA0BBXr3JVVf0Bf6YTsURj1Chcv45Ro3DsGLtbCKtcbItQtuno4NgxTJmCe/eYjiKzYmPh7Ix27XDwIOrerpe13/r1+PABy5YxnYMlxdhCKPOaN8e6dRgwAB/YHdmqLjwc3bphyRIEBbEthtpJOHBm82acOcN0FJa0YgthbTBkCPr2xbBhqMR+I6xPiLB8OSZORHg4hg1jOg1LnIyMcPAgxoxBQgLTUVhSiS2EtcTKlQCwaBHTOWSEcGjM8eOIiUGrVkynYYlfu3ZYsADe3uzAGVYZ2EJYS8jJYd8+7NmDY8eYjiL1kpPh6oqiIly+DCMjptOwJGXiRLi4YPRodsUZVmlsIaw99PRw9CgmTmQnTn3Pgwdo3x49e2L/fqioMJ2GJVkbN+L9e6xYwXQOlpRhC2Gt4uCAP/5Av35s/0/Zjh5F9+5YvRqBgezQmLpIOHBm/XqcPct0FJY0YQthbTNyJFxd2RWHSxMOjZk+HWfOoF8/ptOwmGNsjIMH8cMPeP6c6SgsqcEWwlpo/Xqkp38aPsMCUFgIX1+cOIGYGDg5MZ2GxbT27TFvHvr2zV27dtu2bTvy8/OZTsRiGFsIayEFBRw+jHXrcO4c01GkQFISXF3B4+HSJXZoDOuTiRP5iYk9p0/Pnjz5Q5s2vYjtP6nb2EJYOwknTv3wA54/p/T0dKbjMOb+fbRvDw8P7NvHDo1h/ev58+eKivUFghnFxTNTU/Vfv37NdCIWk9i1Rmut9u0xduw/trZDtbXraWpmRkefMKobDaKkpKQdO/YZGOhqaAz/6SelkBD07ct0JpaUMTQ0JHoG5ACCDx+eX7igP24cO36q7mJbhLVZbOyy4uK1aWnnnj+fvXDhGqbjSEJGRkbr1p7z5xtOmZI8ZsyIs2fZKsgqg5aW1rp1cy0tu1lZ9VyxYtGuXWqtW+PmTaZjsRjCtghrs7y8AkAbAJF2Tk4B03Ek4ebNm1lZvYlGlpRAS6udvT2PvclZZRo8uN/gwZ8GEM+YgT170L8/3N2xahX09ZmNxpI0tkVYmy1aNNXAwFdH5yc5uekODuOZjiN2RUW4cME6P/8mUAC8UFPjy8uzVZBVMQ4Hvr6Ij4eJCZo3R3Awu2xv3cIWwtrMxaV9XNy5U6cGXr78v6Ag25gYpgOJ0+XLcHREQkLjhQuHW1p2trcff/ToZqZDsWSJlhaCghAVhbNn0aoVrl5lOhBLUtjfl2u5evXqdejQAcCWLejfHzduwNyc6UyilpSE2bNx/TqCg9G7N4DR8+aNZjoUS1Y1boxz5xAWhpEj4eyMVatq4bcMqxS2RVhX9OuHqVPRpw/y8piOIjolJQgOhr09TEzw6JGwCrJYIuDlhbg42NrCwQGBgSgqYjrQZxkZGS4u3sbGrdq39/748SPTcWoJthDWIQEBcHbG4MG15PlHVBQcHBAWhqtXERQEZWWmA7FqF1VVBAbi5k3cvQt7e2lZnnTmzOW3bg1JSbkTEzP8l1+CmI5TS7Bdo3XL+vXo2RNz5mD5cqaj1EByMgICcPkyliyBry/TaVi1WoMGCAtDRASmTsXGjVi3DpaWNT3nnTt3liwJMTTUXbToZwMDg28PyMvDmzdITMS7d3j9+tMHwldyc1P4/BEA+Pzmu3efPnkS5uaoXx/168Pc/N+PTUxQ5kCxnJyc4OAtyckf/P1HN2rUqKafSW3BFsK6Rbj6Wrt2sLbGeBkcRsrjYcMGLF6M4cMRHw91daYDseoGd3fExmLjRrRpg4kTMWNG8ezZPyckvF20aHbbtm2rdKrk5OTevSelpa3jct9GRQ3bujUiORkvXiApCV8+yMiAiQmMjWFiAmtrNGuG7t1hbQ1jY8TGjhgxYkJm5hBt7YOhofPbtPnPG6Oi8OIFXrzAmzfQ0IC19ad3Cc9jbIxp04Y9etSrpKTJsWNDY2PP6bMzRQCwhbAO0tXF2bPo0AENGsDNjek0VREdjcmTUa8erlyBrS3TaVh1jIIC/P0xcCB+/RW6uu35/I5Av8jIwXfvnmzZsmVuLrKz//2TkfGfv2ZnIysLmZnIyUFS0v3373sCbQSCNv/8s3r+/BJzc4X69WFri+7dUb8+TE2hq1tuDBMT9+ho01u3brVps6Vp06YAdHRgZ1f6MB4PSUl48wZv3uDtWyQm4to1vH7Nf/TovUAwEUBOzp1bt255enqK8UsmO9hCWBdZW+PAAQwejKgoyETvSEoKZs3CpUtYupTtC2UxydQUu3dj794PwJ8ABIKPbdveLylpqaICTc1Pf7S0oK39719NTNCkCbS1oaUFTU0UFDT38Vn04cNwDudNo0YKly4pVDVD06ZNhSXwO+TlP/WR/pdcgwbc58/vAfULCi7v2TPC2RmGhlW9fi3EFsI6qlMnLF0KDw/cvCml62jk5+erq6sLBAgNRUAABg1CfDw0NJiOxWIB8vIoLn4C2HA4F5cvHzhlCuTkKv9us8OHgxYu/M3QUO+PP0LFF7JM4eF/+fnN/fDh4y+/+MfF2TRrhlmz4O8PRUUJB5EuHJnefyQ0NPTs2bN79+5lOoismjkTt27h4kUoKZXxr4WFhfLy8pJfnGXx4hWBgWsBlXr1jC0solRVORs2lNH5wxItHo/H4/GU2dG3lXDixMlhw2YUFxd37twuIuIQ03GqLyEBc+bg3j0sXQofn0q9JS8vT1VVlSMjK5QLBAI+n6+gUEGzmy2EdZpAgAEDoK2Nv/4q41+ZKoQKCvV5vAeADjBx6FCHffv8JBygbmILYZUQUV5ennqtGK916RJ++gn6+vjzT9jbV3BwrSyE7DzCOo3Lxb59iIvDihVMRwEApKVh3TrweBxACwCHY2xk9ITpUCxWLde1K/7+G4MHo3t3+PoiLY3pQBLHFsK6TkUFJ05g/XocPMhYhoICHD4MLy80bIjoaDRpYs/lunM40+Xlt8ycOZOxWCxWnSEvDz+/T8uO29tj+XIUFzOdSYIYK4RPnz598qTcX/bz8vIePXr0+PHj4jr1v8EQExOcPIkpUyDhVbn5fFy9ivHjYWyMLVvg44N373DoEOLjww4cmDhvnkZKSqyxsbFEM7FYdZiOzqdlx6Oj0bw5wsOZDiQpDIwazcvL692797t37zgcjrGxcXh4eKl+9suXL/fr18/KyqqoqCgjI2Pv3r1du3aVfM46xcFBoqtyx8Vhzx7s2gVdXfj64unT0mO4Bw4c6OHhUTsewLBYsqVxY5w+jYgI+PtjzRqsWVP7h6ox0CIMCQkRCATx8fHx8fFcLnfz5tJ75bRu3TotLe3evXuPHz+ePn26v7+/5EPWQRJYlfv1ayxfjsaNP62OfeUK4uIQEMDOZGKxpI67O+7fh6cnOneGvz8yM5kOJE4MFMKDBw/+8MMP8vLycnJyY8aMOfjNsyk1NTXFz7NarK2t+bVjiWhZIMJVudeu3eLmNnThwpUlJSUfP2LLFnToACcnvHiB7dvx8iWCgtCwoShCs1gs8RAupiN8hNW4MYKD0b27l4aGnbKy5dmqr0H+9u3byMhIHo8n+qA1xkDX6OvXr62trYUfW1tbv3nz5ttjioqKFixYkJWVdfv27Q0bNnznbOnp6REREcKPVVRUXFxcRB64ThHJqtx//bX3t99uZGcvunp1V2joirS0uT16ICAAPXuiomHMLBZLuujpITgYEybA23tDQsIr4HRxcVzv3uMEgsTKn2TKlFmbNh0HrJWUxiYm/q37nUXkmCCWQnjhwoXQ0NIrJsjJyf31118A8vPzv8xVUlZWzs3NLfMkOjo68vLyWVlZsbGxXbp0Ke9az549W7p0qfBjIyMj+wpnwbAq8tdfHDc3VROTYnv7G4aGhuYVPTPk85GczHnzhvvmDffVK87bt9xz565nZ08AGhYX+5eU+P7zT66KCgEoLERhYaUyEFF+fr5MT3KVOcJ5hCUlJUwHkQ117RY1M4OS0iZgNNAMaEa0kssFhwMuFwAUFMDhEIfzaYUaJSXIyYHDgZoaAVBWptjYM0T3AbWCgt87d141YsRvZmZkaSmwshJUOBJg5szf9u27aGysc/z4lgp/HJWSm5urpKRUYd0VSyG0sLDo2bNnqRe53E/dsIaGhl/2k/z48aORkdG3Z1BSUgoICAAwaNCgVq1ajRkzRqOcxbXat2/PTqgXLQ0NnD8PG5uWAIACQ0PtlJQHwn/KyPi0tv2XP0lJePUKSkqf1rm3tkaHDtDVddm2bVNeno6q6q4hQzoYGFR5zAsRcblcdrCMJLET6qukrt2iWVlQUxsO/AV0B+I4nPfPnyMnBxkZAJCSgpISTmEhhD/aMzJQVAQeDxkZHAC5uYiNFXw+E/f5c4OlS5ULC1FSAoEAHA7k5KCkBBUVqKtDSwu6utDXh6EhzMzw8uW+LVv+Jjqdk3PV1XXw+/cPK595ypSATZv2r1wZOmNGp+8fKZZC2Lhx48aNG5f3r61atbp+/bqHhweAa9eutWrV6jun0tTU5PF4AoHgO8ewRC40dBFgDFwG+Kmp9s2b52VkqL1/Dx0dWFnBxgaWlnB2xsCBsLSEhUXphQpzc73Dw4PfvOmtqak6aVIYQ58Ei8USjbAwTJ4MD49ftbVvXbzoraiI06e3W1lV4QyTJvUKCWkJNFJUjHv79u+vW2jv3iEhAS9eIDERSUlITUV6Oh4/xo0byM1FdvZborGAFWCVnh4sbIYKm5vy8uBwoKQEDgcqKgCgrg4OB+rqkJODqiouXz5B9ITPL2sByf9i4Bnh1KlTe/XqZWdnx+FwgoODw8I+/aBs1KjR9u3bO3bsuGPHjvz8/AYNGmRmZq5evXrAgAFaWlqSz1mXpaWlAcJbVQ5Q5PNz9PTU1NWRm4vHj3HzJrS0oKHx6Y+mJnR0/v2rhgaiota9ejWGxxufmnph+vTfjx7dwvDnw2KxqiU1FTNn4sYN7NqFLl0AnKjeEmsbN/4REDA1Li6ue/fupVZtNDWFqSlcXct+47595iNGrCZyAK7o6wtu3UJWFrKykJeHjx9RXIz0dBAhNRXAp48zM8HnIycHRFxAuaCg4ngMFML27dsfOHBg+/btRLRv374OHToIX/f29hZu1mxnZ7dr166zZ89qaWn98MMPY8aMkXzIOm706D82bLACfICPCgoZjx+X7r4W7qz25Y9w9zXhx+np+OefVB6vAwAiu4sX1y9Zgp494eAALruQEYslOw4fxtSpGD0aDx+i5l3mFhYWFhYWVX3XsGHDnj9/tW7daGtrw5Mnz1ZpqlXnzk2jo9upq4cALb9/JLvoNqu0hAR06YJly3iXLvkZGRl9GYtUeTExt3v3npSVNUJD48Tkyf7FxX0jIvDmDTp3hrs7PD1hYlLBGWrTisaygn1GWCW1+xZ9+RLjx+P9e2zbBien//yTbC26ffHixcaNG9f/ZmPGUtj9CFn/kZiI7t0xfz5GjJAfOHBj9baeaN3a+fbtw9HR0Y6O6+0+L0qRkoILF3D6NAICYGICLy+4u8PVlZ1QwWJJEYEA27Zh3jz8/DN++aVK+yxKIzc3t8rMRGcLIetf79+jWzdMmoRx42p6KktLS0tLy69fMTKCry98fcHn4/59hIVh9mw8ewZXV3h5wcPjP0u7vX37NiEhgV1aj8WSpIcP8eOPUFJCdDTKH+9YC7EPbVifZGaiRw8MH45ffhHvheTk4OSEwEDcuYPnz+Hri7t30aYNbGzg74+ICIwaNc3Kyq1bt4Wamg3z8/PFm4bFYgElJVi+HO7uGDsWUVF1qwqCLYQsobw8eHqiUyfMmyfR6+rrw8cHISF4+xb79kFXF7/9ht27TwgEDwWCqJwcz2XLlkk0EItV91y/jpYtcfUq7t6Fnx9k5PGfKLGFkIWCAnh6onFj/PknYxnk5NCmDRYswI0bwsGlwjtTNSWFfYTIYolLfj5mz4aPDwICEBYGMzOmAzGELYR1XUkJBg2CmRm2bpWW3wR9fDy43JZcbm9FxdOnTs3y90c5y/CxWKwqi4qKOn78eH5+/rlzsLXFixeIjYWvL9OxGMUWwjqNz4evL+Tk8NdfUjTJ78CBkEePjh48+ENBwb3YWOWMDLRogaovds9isUobMWJqv367Rox4YGzcddKk/K1bcegQ9PSYjsU0dtRo3UWEiRORno6wMFRrloQYNWnSxNzcnMvlGhpi927873+YOBEbN2LjRknsG8xi1Up8Pv/ChZsZGbcBKCkVr1kT3a1bD6ZDSQWpaQWwJG7mTDx6hOPHRbBmhLh17ox79+DkhFatEBwsgu0SWay6hggnT8plZPCALIDU1f8xNa3HdChpwRbCOmruXEREIDwcsrIyhrIyAgNx7RpOn4azM2JimA7EYsmOGzfg6or58zFjxjJT0y76+i1HjLB1KrVmTB0mZT1iLIlYswZHjuDKFejoMB2liho0wMWLOHwY3t7o0wcrV0JTk+lMLJYUe/wYgYG4dQtz52LsWMjJ9Vy+vPQeeSy2RVjnbNyI9etx+TKqtHytVPHxQXw8lJVha4vdu5lOw2JJpbdvMX48unSBkxOePYOfn8yvlyY+bCGsW0JDERSEixcrXvZayunoIDgYBw5gxQp4euLVK6YDsVhS48MHzJ4NBwfo6ODZMwQEQKniLfnqNLYQ1iEnTmDWLJw7hyrtqCnNOnTAvXvo1g2tWyMwEMXFTAdisRiVl4fly9GkCTIy8OgRgoLA7uVaGWwhrCsiIjB+PE6fhq0t01FESkEB/v64dQsxMWjVCtevo6Cg4NSpU9HR0UxHY7Ekp6QEW7agUSPcvYsbNxASAqPSu4iyysUOlqnljh07GRV1x8qq07Jl3Y4cgaMj04HEw8oKZ84gLAxDhhRkZnYlcldUfOvpeWTXrmCmo7FY4kWEI0cwdy4sLBAWVmu/x8WKLYS12bp123777VJW1kgOZ828eYUdO3oxnUi8vLzA518dPLhTcfFiAOfOtebxeNXbUpHFklrp6elxcXF2dnb16tWLiMDs2eBwsGkT3NyYTiaz2J8Rtdn+/eFZWesAMyLDR49CgFpeCAHUr6+nofHiwwcCcvLzizkc9g5n1Sp37tz19PQrLu7M4cxwdAxJSGi1cCFGjpSWhYJlFPuMsDYzMGjK4ZwE+Coqp1q1asp0HElwdHQcNqyJvn5LIyPXBg2WubsjLY3pTCyW6CxatCk1dXNGxqqPH7elpGx49gy+vmwVrCm2ENZax47h5s25bm6PrazajxiRN3PmZKYTScjatYvT0h4kJ9+7c8fD1RXOzrh1i+lMLFbN8Hi4eBHjxuHCBVXgIwAOJ715czUFdpsyUWA7jmqn9esRFITz59VatNjAdBbGyMkhMBBOTujTB4sXY9w4pgOxWFUkEOD6dRw+jIMHoacHHx9ERs4aNWpAbu5KNbWc5cuPMB2wlmALYW1DhIULcfAgrl9H/fpMp5ECXl64cgUDBiA6GiEhUFFhOhCLVZEv9e/QIejqwscHV6+iQQPhP5olJNzKysrSYmcIio5sd43eunUrhl19+St8PiZMwJkzuHJFBFVwzpw5O3bsEEWuKrt48eKoUaNEdbZGjXDjBoqK0KEDXr4U1Vlrm7/++uvXX39lOoXMiIyMHDlypGjPKRDg6lX4+8PMDOPHQ0cHV64gLg6BgV+q4CcMVsHRo0efP3+eqatX1dmzZyvzk0S2W4Q8Ho/H4zGdQloUFWHECGRkIDISGhoiOGFhYWFhYaEITlR1JSUleXl5IjyhujoOHsSWLWjfHrt2oXt3EZ67lmDwv1sW1eQWLSwsjIiI0NbW7tChA75q/x0+DB0d+Pjgf/9Do0YijSs6eXl5JSUlTKeorOLi4sr8N8l2IWR9kZkJb2+Ym+PsWbDPz8vj54emTTF0KH78EfPngyvbHSIsmVRQUODg4J6U5Kqg8P/27j2oiWsNAPjJg0cgGPKAIEQjUoRiFJmq0GJ4BV8IBtRppdJSVIrSjo5MO2CZzjjSitVSadXqRS2Eqx1LLxT7UkqJiYGolWujQQQvEFQqKhCCCRiEZO8fO83NRUBEJGC+3185jz35jGS/7O7ZPX8tWFDk4/PVv/6FnJ0nev57sUEifBHcvYuWL0d8PsrNhZ37E/D56OJF9Prr6N//Rv/8J3J2tnRAwDr09qLbt9GtW+j06aqbN4P1+l0IoYqKwEWL+qVS8oAzn2CcETAMs3QMo5eVlXXgwIG5c+daOhALq6vbSqG0crlFYzusSqWiUChulnhkoUajuXv3rq+v73Ma32gkNzS8SyAYvL3/8ZzeYtK5d+9eT0+P5wvzRPbn6eHDh5cv3+7tRV5eDlwuZ5iefX3UpqZ39Hp2by8DIQKJ1E0iXXzw4D9E4odkcjeB8Hpg4CRbFLS+vp7NZjtPkp+QnZ2d/v7+x44dG77b5E6ECCG5XN7T02PpKAAAAExEgYGBTk+aNDHpEyEAAADwLOCCEgAAAKsGiRAAAIBVg0QIAADAqkEiBAAAYNVIO3bssHQMo9HU1HT16tWbf2tvb3d3d7d0UC+Itra2S5cuzZgxw1Qjk8lsbW2fOPNqJO7cuXP9+nUPDw9TjUQicXZ2tre3N+926tSpK1eu8Hg880qDwSAWi93d3fG1dvV6vUQisbe3H5PArNmNGzdqampM3yaNRmORe2Ymi9u3b9++fdvV1dXSgYyTmpqaBw8eMBgMvCiTyTAMM90+cfbsWTqdPuD7a3H4voLJZJoCu3btWltb25D/a9jklJaWNnXq1Ii/bdiwwdIRvTj0er2fn19BQQFeLC0t9fDw6OzsHJPBjx07tnDhQvOaKVOmVFZWDuiWnZ2dmZk5oFKr1SKEmpubMQzT6XSLFy+OiYl5+PDhmARmzTZs2DBt2jTTt2nLli2WjmhCy87Ojo6OtnQU4ycrK2vlypX4a41GQyaT33rrLbzY0tJCIBDu379vueiGlJSUlJCQgL9uaWlhMpmP72dMJvGTZZYuXZqfn2/pKF5AdnZ2x44dW7FiRUREBIVCSUlJEYlEE+r+2c7OzhUrVsycObOgoAA/OgTPaNWqVbm5uZaOAkxEYWFhOTk5BoOBRCLJZLIlS5bIZDK8SSKR8Hg8FxcXy0Y4qH379s2ZM6ekpCQuLm7jxo3JycnBwcFDdYadCBhEUFBQUlLSpk2bHBwcYmNjly5daumI/ufevXtCoTA4OHj//v1EeKAcAM/ZwoUL+/r6lErlvHnzpFJpTExMa2trU1PTzJkzpVJpWFiYpQMcHI1Gy8vLW79+/Y0bN27dulVaWjpMZ9iPgMFlZWVdu3btjz/+2Lt3r6Vj+T8xMTHh4eEHDx6ELAjAOLC1tX311VclEglCSCKRhIaGhoSESKVShNBEToQIoWXLloWHh2dmZubn59vZ2Q3TE3YlYHB1dXVdXV06na6rq8vSsfyfBQsWlJWV3bt3z9KBAGAtQkNDpVKpVqv966+/fH198UTY2tra0NAQEhJi6eiG1N3dffHiRUdHx5aWluF7QiIEg+jv79+4ceOuXbuSkpJSUlLGcGQnJyd8zgvOYDD09PQ81bTPAwcOhIWFhYaGtra2jmFgAIChhIWFnTt3TiqVBgcHEwiEkJCQs2fPSqVSHo/HYrEsHd2QMjIyZs+eXVJSsnnz5ra2tmF6QiIEg9i5cyeNRtu0adOnn36qUqkKCwvHamQ/P7+mpibT8dylS5dsbGxeeppFaAgEwsGDByMiIsLDwyEXAjAOFi5c2Nvbe+DAgdDQUIQQi8WiUqkFBQUT+bxoVVVVUVFRXl5eZGTkihUr0tLShukMk2XAQAqFYv/+/dXV1QQCAZ9BGh0dLRAIzG/+G7XZs2fHx8cvXrz47bfffvTo0aFDh3bs2OHg4PBUg+C5MDU1NSIiQiwWT5069dkDA2DklErl5s2b8dc0Gm337t2Wjed5wy8TlpWV7dmzB68JCQk5fPjw2J4uGkPd3d3vvPPOl19+yWazEUL79u3j8XilpaWxsbGD9p+sN9Q7OzvPmTOHy+VaOpAXUPiESBEAAAhrSURBVG1t7ZtvvhkQEIAXp02b9vLLL5NIpLHKN7GxsVwu9+bNmyQS6YMPPli7du3jfSorKx89ehQREWFeSSAQOBxOYGCgra0tgUCIioqi0+n9/f3m9/6DUaDT6f7+/hzOcOvqARMnJycvL6+pf+NwONawJKq3t3dgYOCyZcsIBAJCyNPT09/fPyYmZqLdSo9TqVQ+Pj7x8fF40d7ePjQ0tLu728fHZ9D+sAwTmIh2796t0+k++eQTSwcCAHjxwTVCAAAAVg2OCMFEdP/+fYPBABf/AADjABIhAAAAqwanRgEAAFg1SIQAAACsGiRCAAAAVg0SIQAAAKsGiRAA8AQlJSXl5eWWjgKA5wVmjQIAnmD+/Pmenp7ff/+9pQMB4LmAI0IAxk97e7terx/dthqNZqwWn2pra1Or1UO1qtXqjo6OEQ7V2dl5//794Tt0dnY+XXwAjC9IhACM0unTpxkMxo0bN/Di0aNHGQyG6eG9Wq3WxcXlyJEjCKGOjo64uDgnJycXFxdHR0dfX9/i4mK8W1NTE4PByM/PNx+5pKSEwWD8+eefePHUqVN+fn50Ot3NzW369OknT54cNJ7ffvuNwWCcO3fOvHLPnj1sNtuU9g4fPjx9+nRXV1cmk8nj8fDVVk2OHDkyc+ZMJpPJYrGYTObnn3+OEJo3b55Cofjxxx8ZDAaDwUhISMA7l5WVzZ49m8FgsNnsGTNmfPvtt6ZxcnNzGQzGpUuX/P39GQzG6tWrn+6TBWCcYQCAUVGr1SQS6euvv8aLa9assbW1XbBgAV78+eefEUJKpRLDsObm5uTk5F9//bW2traysnLVqlVkMvny5ct4z6CgoODgYPORo6KiZs2ahb/+4YcfiERiYmLi+fPnL1++nJqaSiAQzpw583g8fX19bm5uSUlJ5pU+Pj5CoRB/vWfPHiKRuH379urq6gsXLgiFQgqFUltbi7fu3bsXIbRu3TqJRKJUKk+cOLF3714Mw2Qymbe3N5/PLy8vLy8vVygUGIZdvHjRxsZm0aJFEolELpfHxsYSCITi4mJ8qF27diGEPD09c3Jy5HK5WCx+xo8agOcKEiEAozd//vzVq1djGGYwGFgsVmpqKolEUqvVGIZt27aNzWYbjcbHt+rr6+NwOOnp6Xjx0KFDCKG6ujq8ePfuXTKZvHv3bgzDjEajl5fXkiVLzDcPCQmJjIwcNJ5t27ZRqVStVosXq6qqEEIlJSUYhmk0GkdHx/fff9/UWa/Xc7nczZs3YxjW1dVFpVKjo6MHHfaVV15Zs2aNeY1QKHR2du7q6sKL/f393t7ec+fOxYt4Ijx06NBQnxsAEwqcGgVg9AQCgVgsNhgMCoWio6Nj+/btjo6O+PnGiooKgUCAr1mDENJqtUeOHElPT09JSXnvvfeMRmNjYyPeFB8fT6FQTpw4gRePHz9uNBrXrVuHEGpqampsbJw1a9bvZjgcjlKpHDSe9evX63S60tJSvCgSiZhMZlRUFEKoqqqqu7ubw+GYxpHJZFwut6amBiF04cIFnU63YcOGEf7DFQrF8uXLp0yZghdJJNKaNWuUSqX55UChUDjyTxIAC4KFeQEYPYFA8NlnnykUCrFYzOPxOBwOn8+vqKjg8/k1NTVbt27Fu125ckUgEFAolMjISCaTSSaTbWxsurq68FYajSYUCgsLC3fs2EEkEgsLC5cuXYqvDojPjikoKDClSROj0UgkDvwhy+PxAgICRCJRQkKCXq8vKipKTEy0s7MzDZWdnT1gK09PT4RQe3s7QmiESxIajcaWlpYBj0R3d3fHMEytVtPpdLwGXxMVgIkPEiEAo8fn8+3t7SsqKsRisUAgQAgJBIK8vDw+n280Gk0LC+fm5k6ZMkWpVDo6OuI1v/zyi/k4iYmJJ0+elEqlNBrt6tWrH330EV5Po9EQQjk5Oe++++4IQ0pMTExLS7t165ZcLtdoNImJieZDlZaWhoWFPb6Vs7Mz+jtZPhGRSKRQKHjuNGlrazO9i6nbCGMGwLLgLxWA0bO3tw8KCjp9+nRlZaUpEdbV1YlEIi8vrxkzZuDdVCqVr6+vKQs2NDTU19ebj7NkyZJp06aJRCKRSESj0VauXInX+/r6urq6PtUNfAkJCTY2NsePHxeJRPgBIl7/2muv2djYFBUVDbpVYGCgnZ3dUK1UKvXhw4fmNUFBQb///ntvb6+p5qeffnrppZdYLNbIQwVgorD0RUoAJresrCyEEJlMxmeOGI1GV1dXhFBycrKpz9atWykUypkzZ/R6fXV1dUBAAJVKXbx4sfk4GRkZVCqVxWJt2rTJvD4vLw8htH79+rq6up6ensbGRpFIhE+lGUpcXByHwyGRSF988YV5/YcffkgkEj/++GOVStXT01NXV/fVV1/l5+fjrenp6UQiMTMzU6VS6XQ6uVwuEonwptTUVDqdXlxcXF1d3dDQgGFYWVkZgUB44403mpub79y5s2XLFoRQXl4e3h+fLDOqjxMAC4A/VgCeiVwuRwiZ3/+wdu1ahNB3331nqmlvbw8ODsZ/etrZ2WVlZUVGRg5IhKZjxPPnzw94i6NHj7q5uZl+vLq4uAzIcAPgk2XIZHJra6t5vcFg2Llzp2mGC0KIy+UWFRXhrf39/ZmZmRQKBW8ik8kZGRl4U0tLS1RUFH76NDY2Fq/85ptvTJcDHRwcsrOzTW8EiRBMLvCINQDGA4ZhKpVKrVbPmjXLPBWNkNForK+v12q1bDYbP9obdSR9fX3Xr1/v7e11d3f38PAY0IofKRKJRC6Xa8pzQ3n06FFtbW1/f7+fn5+Dg8OoQwLAsiARAgAAsGowWQYAAIBVg0QIAADAqkEiBAAAYNUgEQIAALBqkAgBAABYNUiEAAAArNp/ASZAb/RqfEZfAAAAAElFTkSuQmCC", "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 [-1.0321137134594964e-16, -6.882802718489987e-16, -4.437681174180191e-16]\n [-8.674337307691979e-16, 4.962781517099956e-16, -1.0005301709361671e-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 almost zero in this highly symmetric configuration." ], "metadata": {} }, { "cell_type": "markdown", "source": [ "## Where to go from here\n", "Take a look at the\n", "[example index](https://docs.dftk.org/stable/#example-index-1)\n", "to continue exploring DFTK." ], "metadata": {} } ], "nbformat_minor": 3, "metadata": { "language_info": { "file_extension": ".jl", "mimetype": "application/julia", "name": "julia", "version": "1.7.3" }, "kernelspec": { "name": "julia-1.7", "display_name": "Julia 1.7.3", "language": "julia" } }, "nbformat": 4 }