{ "cells": [ { "cell_type": "markdown", "id": "edf16ac1-0706-4779-b576-6abefc4d28d4", "metadata": {}, "source": [ "# Nuclear Quantum Effects With RPMD and adQTB\n", "\n", "Nuclear quantum effects (NQEs) play an important role in many molecular systems. This is especially true for low temperature simulations, but in some cases they can even be significant at room temperature. For example, they are important in systems like water where hydrogen bonding plays a major role (https://doi.org/10.1021/acs.chemrev.5b00674).\n", "\n", "OpenMM has two methods for simulating nuclear quantum effects: Ring Polymer Molecular Dynamics (RPMD) and Adaptive Quantum Thermal Bath (adQTB). This tutorial shows how to use both of them and compares the results.\n", "\n", "## Example: Parahydrogen\n", "\n", "Molecular hydrogen ($H_2$) can exist in two nuclear spin states. In orthohydrogen the nuclear spins are parallel to each other, and in parahydrogen they are antiparallel. Parahydrogen is the more common state at low temperature, while orthohydrogen is the more common state at high temperature.\n", "\n", "\n", "\n", "For this tutorial, we will simulate a box of parahydrogen molecules at 25 Kelvin. First let's import the libraries we will be using." ] }, { "cell_type": "code", "execution_count": 1, "id": "bedbc832-2673-4ac5-863e-fe99dc2ee9ab", "metadata": {}, "outputs": [], "source": [ "from openmm import *\n", "from openmm.app import *\n", "from openmm.unit import *\n", "import numpy as np\n", "import matplotlib.pyplot as plot" ] }, { "cell_type": "markdown", "id": "b9864577-a8a6-4063-9c62-e107d7d38172", "metadata": {}, "source": [ "Now we build a system to represent the box of molecules. Each molecule is represented by a single particle. They interact with each other through an empirical pairwise potential (https://doi.org/10.1063/1.437103)." ] }, { "cell_type": "code", "execution_count": 2, "id": "786404c2-6e74-491c-a7ee-578c3f8304f8", "metadata": {}, "outputs": [], "source": [ "particles = 32\n", "box_size = 1.1896\n", "temperature = 25*kelvin\n", "system = System()\n", "system.setDefaultPeriodicBoxVectors(Vec3(box_size, 0, 0), Vec3(0, box_size, 0), Vec3(0, 0, box_size))\n", "force = CustomNonbondedForce(\n", " \"\"\"2625.49963*(exp(1.713-1.5671*p-0.00993*p*p) -\n", " (12.14/p^6+215.2/p^8-143.1/p^9+4813.9/p^10)*(step(rc-p)*exp(-(rc/p-1)^2)+1-step(rc-p)));\n", " p=r/0.05291772108; rc=8.32\"\"\")\n", "force.setNonbondedMethod(CustomNonbondedForce.CutoffPeriodic)\n", "force.setCutoffDistance(box_size/2)\n", "system.addForce(force)\n", "for i in range(particles):\n", " system.addParticle(2.0*amu)\n", " force.addParticle()\n", "positions = np.random.rand(particles, 3)*box_size" ] }, { "cell_type": "markdown", "id": "80cc95d1-94e1-46ca-934c-da7ca5b795e0", "metadata": {}, "source": [ "We are going to simulate this box of parahydrogen and record the pairwise radial distribution function (RDF). The following function runs the simulation and computes the RDF." ] }, { "cell_type": "code", "execution_count": 3, "id": "0fb7b29d-0d1e-4afe-89d8-0e457097d2c1", "metadata": {}, "outputs": [], "source": [ "def compute_rdf(context):\n", " bins = 100\n", " iterations = 2000\n", " counts = [0]*bins\n", " for _ in range(iterations):\n", " # Run a few steps of dynamics.\n", "\n", " context.getIntegrator().step(20)\n", "\n", " # Count the number of distances in each bin.\n", "\n", " pos = context.getState(positions=True).getPositions().value_in_unit(nanometer)\n", " for i in range(particles):\n", " for j in range(i):\n", " delta = pos[i]-pos[j]\n", " delta -= np.round(delta/box_size)*box_size # Apply periodic boundary conditions\n", " dist = norm(delta)\n", " counts[int(bins*dist/box_size)] += 1\n", "\n", " # Convert the histogram of distances to the RDF.\n", "\n", " scale = (box_size*box_size*box_size)/(iterations*0.5*particles*particles);\n", " rdf = []\n", " for i in range(bins//2):\n", " r1 = i*box_size/bins\n", " r2 = (i+1)*box_size/bins\n", " volume = (4.0/3.0)*np.pi*(r2*r2*r2-r1*r1*r1)\n", " rdf.append(scale*counts[i]/volume)\n", " return rdf" ] }, { "cell_type": "markdown", "id": "56b1b498-e1ad-4b83-927a-204218a0ed67", "metadata": {}, "source": [ "Before we look at NQEs, first let's run a simulation with a standard Langevin integrator. This gives us the classical distribution." ] }, { "cell_type": "code", "execution_count": 4, "id": "e52c2f16-d5ad-4b8c-baaa-bd47206d9387", "metadata": {}, "outputs": [], "source": [ "integrator = LangevinIntegrator(temperature, 1.0/picosecond, 1*femtosecond)\n", "context = Context(system, integrator)\n", "context.setPositions(positions)\n", "LocalEnergyMinimizer.minimize(context)\n", "context.setVelocitiesToTemperature(temperature)\n", "integrator.step(1000) # Equilibrate before collecting data\n", "classical_rdf = compute_rdf(context)" ] }, { "cell_type": "markdown", "id": "2ae983b0-7990-4671-bec4-4f3a820556d8", "metadata": {}, "source": [ "Here is a graph of the RDF, so we can see what it looks like." ] }, { "cell_type": "code", "execution_count": 5, "id": "dbfc0367-04a3-4071-9064-6badfd05b885", "metadata": {}, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAiMAAAGdCAYAAADAAnMpAAAAOnRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjEwLjMsIGh0dHBzOi8vbWF0cGxvdGxpYi5vcmcvZiW1igAAAAlwSFlzAAAPYQAAD2EBqD+naQAAN3VJREFUeJzt3Xt8lPWd9//3ZCZnkpADOUHAcApnxIASKoiiUVC3VrvVbbt693DfN/1p1bLetuhubbW9cbv+XGqtcrue6u222l3QuhWVtEJQQQUMcj4HEkImISFkQsh5rvuPZIYEAmTCzFxzzbyej8c8aq7MZD65Csmb7+HztRmGYQgAAMAkUWYXAAAAIhthBAAAmIowAgAATEUYAQAApiKMAAAAUxFGAACAqQgjAADAVIQRAABgKofZBQyE2+3WsWPHlJSUJJvNZnY5AABgAAzDUFNTk3JzcxUVdf7xD0uEkWPHjikvL8/sMgAAwCBUVlZqxIgR5/28JcJIUlKSpO5vJjk52eRqAADAQLhcLuXl5Xl/j5+PJcKIZ2omOTmZMAIAgMVcbIkFC1gBAICpCCMAAMBUhBEAAGAqwggAADAVYQQAAJiKMAIAAExFGAEAAKYijAAAAFMRRgAAgKkIIwAAwFSEEQAAYCrCCAAAMBVhBEBE++vuGv1xU6UMwzC7FCBiWeLUXgAIhIbmdi1+fYs6ugx9UdGgX35tquxRFz5dFID/MTICIGKt2eVUR1f3iMgbmyp1779/obbOLpOrAiIPYQRAxPrztmpJ0rzxwxRjj9L7O5367qubdKqt0+TKgMhCGAEQkU40t2vDwXpJ0s//ZrJe+c4sJcbY9cmBen3z3z7VieZ2kysEIgdhBEBE+mCnU11uQ5Nzk5WfkaivjM3Q7//7bKUmRGvb0Ub97YoNOnayxewygYhAGAEQkVZv756iWTQ1x3ttet5Q/cfiIuWkxOng8WZ9/fkNOnj8lFklAhGDMAIg4tSfavNO0dzcK4xI0tjMJP3nD+Zo9LBEHWts1d+u2KjtRxvNKBOIGIQRABHng5016nIbmjI8WZdlJJ7z+eFD4/Uf/7NIU4en6ERzu+56YaM2HKwzoVIgMhBGAESc/qZozpY+JFZ/+B+zVTQ6Xc3tXfpvr2xiygYIEMIIgIjSPUXTPcpx9hTN2YbEOvTKd2Zp5qhUtXe69W7PVmAA/kUYARBR3t/plNuQpg5P0aj0c6dozhYXbdftV4yQJK3fdzzQ5QERiTACIKJ4pmhunnbhUZHe5o3PkCR9UdGgxtMdAakLiGSEEQARo+5UmzaeZxfNhYxITdCYYYlyG9InLGQF/I4wAiBivL+je4pm2ogU5aUl+PTaa8ZnSmKqBggEwgiAiOGdovFhVMTDM1VTuu+4DMPwa11ApCOMAIgIx5va9Omh7imaC23pPZ/Zo9MV64hSdWOr9teyxRfwJ8IIgIjg2UUzfRBTNFL3rpqrRqdLYqoG8DfCCICIsHqb77tozjZv3JmpGgD+QxgBEPZqm1r1WXn3FM3CKYMPI/MLhkmSPjt0QqfbO/1SGwDCCIAI8EHPLprpeUMHNUXjMWbYEA0fGq/2Lrc+O3TCjxUCkY0wAiDsvduzi+aWQSxc7c1ms2ne+O7REaZqAP8hjAAIa91TNN2jGAunZl/y17umZ4svi1gB//EpjCxbtkyzZs1SUlKSMjMzddttt2nv3r0XfV1paakKCwsVFxen0aNHa8WKFYMuGAB88f4OpwxDujxvqEakDn6KxmPO2AzZo2w6VNesyhOn/VAhAJ/CSGlpqe699159+umnKikpUWdnp4qLi9Xc3Hze15SXl2vRokWaO3euysrK9Mgjj+j+++/XypUrL7l4ALiYP/fsornlEnbR9JYcF63CkamSmKoB/MXhy5Pff//9Ph+/8soryszM1JYtWzRv3rx+X7NixQqNHDlSy5cvlyRNnDhRmzdv1lNPPaU77rhjcFUDwADUulq16bBnisY/YUTq7sb6+eETKt13XN+ePcpvXxeIVJe0ZqSxsVGSlJaWdt7nbNy4UcXFxX2u3Xjjjdq8ebM6Ovo//bKtrU0ul6vPAwB89V7PFM2MkUM1fGi8376u55yaDQfq1N7p9tvXBSLVoMOIYRhasmSJrr76ak2ZMuW8z3M6ncrKyupzLSsrS52dnaqr6//0y2XLliklJcX7yMvLG2yZACLYu9sGfxbNhUzOTVZ6Yoya27v0RUWDX782EIkGHUbuu+8+bdu2TX/4wx8u+lybzdbnY88hU2df91i6dKkaGxu9j8rKysGWCSBCHW9q06Yj3VM0gzmL5kKiotjiC/jToMLID3/4Q73zzjtau3atRowYccHnZmdny+l09rlWW1srh8Oh9PT0fl8TGxur5OTkPg8A8MUep0uGIY0ZlqhcP07ReHhP8d1LGAEulU9hxDAM3XfffVq1apU+/PBD5efnX/Q1RUVFKikp6XNtzZo1mjlzpqKjo32rFgAGqKJn2+1l6YkB+fpzx3WPjOyqdqm2qTUg7wFECp/CyL333qvXX39dv//975WUlCSn0ymn06mWlhbvc5YuXaq7777b+/HixYt15MgRLVmyRLt379bLL7+sl156SQ899JD/vgsAOEvlie6fS5fS/v1CMobEaurwFEnSR/v6X/8GYGB8CiPPP/+8GhsbNX/+fOXk5Hgfb775pvc51dXVqqio8H6cn5+v1atXa926dbr88sv1xBNP6JlnnmFbL4CA8jQkC1QYkXpN1bBuBLgkPvUZ8Sw8vZBXX331nGvXXHONvvjiC1/eCgAuiWeaZmQAw8g14zP127UH9dH+4+pyG7JH9b8oH8CFcTYNgLAUjDAyY+RQJcU61HC6QzuqGgP2PkC4I4wACDuNLR1qbOluqjgi1f87aTyi7VGaM7Z7VyBTNcDgEUYAhB3PepGMITFKjPVpNtpnnm6shBFg8AgjAMJOMBavengWsZZVNKjxdP9HXAC4MMIIgLATjPUiHiNSEzQ2c4jchvTJQbb4AoNBGAEQdiobekZGUgMfRiRpXk8DNLqxAoNDGAEQdip6Gp4FY2REkq4p6A4j6/cfH1ALBAB9EUYAhJ1grhmRpKvy0xTriFJ1Y6v2154KynsC4YQwAiCsdLkNVTX0jIykByeMxEXbddXo7i2+6/bWBuU9gXBCGAEQVmpcrWrvcivablN2clzQ3nf++O6pmrV7WDcC+IowAiCseHbSDB8aH9T27Asmdvcb2XT4hFytbPEFfEEYARBWKoK8XsRjVHqiRg9LVKfb4BRfwEeEEQBh5ahJYUSSFkzoHh35656aoL83YGWEEQBhJZgNz852bU8YKd3bfYovgIEhjAAIK2aGkVmXpSkp1qH65nZ9efRk0N8fsCrCCICwEuyGZ71F26M0z7urhi2+wEARRgCEjZb2LtWdapMUvFbwZ7vOs25kN2EEGCjCCICw4TmTJjnOoZSEaFNqmF8wTDabtKvaJWdjqyk1AFZDGAEQNirqe9aLBKnzan/Sh8Tq8ryhkqQPmaoBBoQwAiBseHuMmDRF43FdQfdUDWEEGBjCCICw4ZmmMWPxam+eLb6fHKhTa0eXqbUAVkAYARA2gn1a7/lMzk1WVnKsWjq69OmhelNrAayAMAIgbJjZY6Q3m83m3VXDFl/g4ggjAMKCYRiq7OkxYvbIiCRdNyFLkvTXPbUyDLqxAhdCGAEQFupOtaulo0s2W/eJvWb7yth0xTiidLShRQdqT5ldDhDSCCMAwoJniiY3JV4xDvN/tCXEOFQ0Ol1S9+gIgPMz/28sAPjBmcWr5o+KeHjWjbDFF7gwwgiAsFAZIj1GevOEkS1HGtR4usPkaoDQRRgBEBZCZSdNb3lpCRqXOURdbkOl+4+bXQ4QsggjAMKCN4yY2Aq+P9dN7Jmq2V1jciVA6CKMAAgLnmmaESE0TSOdaQ2/bt9xdbnZ4gv0hzACwPLaO92qdnWfkBtK0zSSVDgqVclxDp083aGyigazywFCEmEEgOVVnWyRYUjx0XZlDIkxu5w+HPYoXcPBecAFEUYAWF7vxas2m83kas513YRhkggjwPkQRgBYXkUI9hjp7ZrxmYqySXucTao62WJ2OUDIIYwAsLyjIXJa7/mkJcZoxshUSYyOAP0hjACwvFDsMXI2TvEFzo8wAsDyKkKw++rZPGHkkwN1amnvMrkaILQQRgBYXqg2POttQnaSclPi1Nbp1sZDdWaXA4QUwggAS2s83aGm1k5JoT0yYrPZdG3P6MhfdzNVA/RGGAFgaZ5RkWFJsYqPsZtczYV5pmrW7T0uw6AbK+BBGAFgaWfWi4Tmtt7eisakK8YepaqTLTpQe8rscoCQQRgBYGlW2EnjkRDj0FWj0yR1j44A6EYYAWBplQ3WCSOSdM347m6s6/axbgTwIIwAsLTKEG94drb5PefUbCpvUHNbp8nVAKGBMALA0iosFkbGDEvUiNR4tXe5teFgvdnlACGBMALAsrrchqoaus96sco0jc1m0/yCnqmavUzVABJhBICFVTe2qNNtKMYepazkOLPLGbBrC9jiC/RGGAFgWZ4pmuGp8bJH2UyuZuB6b/E9eJwtvgBhBIBlWW3xqgdbfIG+CCMALKvyhGe9SOg3PDubd4svYQQgjACwLis1PDubZ4vv5+Un2OKLiEcYAWBZZ1rBWy+M9N7iu5EtvohwhBEAlmXVNSPSWVt86caKCEcYAWBJzW2dqm9ulySNTLdeGJGk+ePZ4gtIhBEAFuU5kyYlPlrJcdEmVzM4c8Z2b/E92tCig8ebzS4HMA1hBIAlVdRbd/GqR0KMQ1fme7b4MlWDyEUYAWBJVt5J09uZ1vBs8UXkIowAsKSjPWfSWHHxam+eMMIWX0QywggASzpzWq/1Gp71NmbYEA0fyhZfRDbCCABLCpdpGpvNpmsnsMUXkY0wAsByDMPw9hixehiR2OILEEYAWM7xpja1dboVZZNyh1p7mkZiiy9AGAFgOZ4pmpyUeEXbrf9jjC2+iHTW/1sMIOJ4Gp5ZffFqb55dNaX72OKLyEMYAWA5ta42Sd0jI+HCE0Y+O3RCp9vZ4ovIQhgBYDmeM2nSE2NMrsR/2OKLSOZzGFm/fr1uvfVW5ebmymaz6e23377g89etWyebzXbOY8+ePYOtGUCEqzvVPTKSPiTW5Er8p88pvnRjRYTxOYw0Nzdr+vTpevbZZ3163d69e1VdXe19jBs3zte3BgBJUv2pnpGRIeEzMiJJ8wt6tvjuq2WLLyKKw9cXLFy4UAsXLvT5jTIzMzV06FCfXwcAZ6tv7h4ZyQizMDJnTPcW38oTLTpU16wxw4aYXRIQFEFbMzJjxgzl5ORowYIFWrt27QWf29bWJpfL1ecBAB7ekZHE8JmmkaTEWIdm5adKYqoGkSXgYSQnJ0cvvPCCVq5cqVWrVqmgoEALFizQ+vXrz/uaZcuWKSUlxfvIy8sLdJkALMIwjDMLWMNsZETq3Y2VfiOIHD5P0/iqoKBABQUF3o+LiopUWVmpp556SvPmzev3NUuXLtWSJUu8H7tcLgIJAEnSqbZOtXe6JYXfyIgkXTshU79cvVufHqpXY0uHUuKjzS4JCDhTtvbOnj1b+/fvP+/nY2NjlZyc3OcBANKZKZqEGLviY+wmV+N/YzOHaHzWEHV0GSrZVWN2OUBQmBJGysrKlJOTY8ZbA7A4z+LVcJyi8bh5aq4k6d1tx0yuBAgOn6dpTp06pQMHDng/Li8v19atW5WWlqaRI0dq6dKlqqqq0muvvSZJWr58uS677DJNnjxZ7e3tev3117Vy5UqtXLnSf98FgIhRF6aLV3u7eVq2/vUv+/TR/jo1nu5QSgJTNQhvPoeRzZs369prr/V+7Fnbcc899+jVV19VdXW1KioqvJ9vb2/XQw89pKqqKsXHx2vy5Ml69913tWjRIj+UDyDSeKZpwm1bb29jM5M0ITtJe5xN+mCXU9+YyZo5hDefw8j8+fMv2Izn1Vdf7fPxww8/rIcfftjnwgCgP/We7qthPDIiSTdPzdEeZ5Pe3VZNGEHY42waAJYSztt6e1s0rXtd3ScH6tTQ8z0D4YowAsBSwvFcmv6MGTZEE3OS1ek2tGaX0+xygIAijACwlEhYM+JxS8/oyLvbCSMIb4QRAJbi3dob5mtGJGnRVKZqEBkIIwAsJVxP7O1PfkaiJucmq8tt6IOdjI4gfBFGAFhGl9tQw+nICSOSdLN3qqba5EqAwCGMALCMk6fb5e7pLJCaECFhpGeqZsPBeu+2ZiDcEEYAWIZnW+/QhGhF2yPjx9eo9ERNGe6ZquGsGoSnyPjbDCAseLf1JkbGqIiH96ya7ZxVg/BEGAFgGWcWr4b/TprePFM1Gw/WewMZEE4IIwAsw7NmIhJ6jPQ2Mj1B00akyG1I7+9gVw3CD2EEgGV4W8FHQI+Rs3lGR97dxq4ahB/CCADLqIugHiNn8zRA+6y8XrVNrSZXA/gXYQSAZdRHyLk0/clLS9D0vKFyG9IHTNUgzBBGAFiGZ5omI8J203jc0jM68memahBmCCMALCOSR0YkaeHUbEnS54dPqNbFVA3CB2EEgGVE0rk0/RmRmqAZI4fKMKT3mKpBGCGMALCEts4uNbV1SpIyInA3jQe7ahCOCCMALOFEz3oRR5RNyfEOk6sxj2dXzaYjJ+RsZKoG4YEwAsASPFM0aYkxstlsJldjntyh8SocldozVcPoCMIDYQSAJdRF+OLV3piqQbghjACwBM/ISKS1gu+PZ6pm85EGVZ1sMbka4NIRRgBYQn1zZJ7Y25/slDjNHp0mSXq7rMrkaoBLRxgBYAmRemLv+dxxxQhJ0sotR2UYhsnVAJeGMALAEiL5XJr+LJyao/houw7VNWtr5UmzywEuCWEEgCV4pmkiucdIb0NiHbppSndH1pVfHDW5GuDSEEYAWEKkd1/tj2eq5r++rFZbZ5fJ1QCDRxgBYAmRfi5Nf4rGpCs7OU6NLR36cHet2eUAg0YYARDyDMPwntjLbpoz7FE2fe2K4ZKYqoG1EUYAhLzm9i61dbolMU1ztjt6wsi6vce9jeEAqyGMAAh5nimahBi7EmIi91ya/ozNTNL0ESnqdBt6Z+sxs8sBBoUwAiDk1fU6lwbnuqOwp+cIUzWwKMIIgJDH4tULu3VarqLtNu085tIep8vscgCfEUYAhDzP4tUMRkb6lZoYo+smZEqSVn1Be3hYD2EEQMg7MzJCGDkfT8+Rt8qq1NnlNrkawDeEEQAhr45zaS5qfkGm0hJjdLypTR8dqDO7HMAnhBEAIY8eIxcX44jS30zPlcRUDayHMAIg5HmmaTIYGbkgz1TNmp1OuVo7TK4GGDjCCICQx7k0AzNleLLGZw1RW6db726rNrscYMAIIwBCnufE3nRO7L0gm82m23tGR1bRcwQWQhgBENLcbkMnPFt7GRm5qK/NGK4om7TpcIOO1DebXQ4wIIQRACHtZEuH3Eb3f6eygPWispLjdPW4YZJYyArrIIwACGmexatDE6IVbedH1kB4Ds9bVXZUbk+SA0IYf7MBhDTOpfFd8aRsDYl1qPJEizYdPmF2OcBFEUYAhDTP4tUMFq8OWHyMXTdPzZHEVA2sgTACIKSxrXdwbu+Zqnl3e7Va2rtMrga4MMIIgJDGuTSDM+uyNOWlxetUW6dWlB40uxzggggjAEJanbcVPNM0voiKsumh4gJJ0m8+3K/NrB1BCCOMAAhpZ1rBMzLiq69ePly3zxgutyE98MZWNbbQIh6hiTACIKTVc2LvJfn5VydrZFqCqk626NG3tssw2OqL0EMYARDSOLH30iTFRevXd10ue5RNf95WrZXsrkEIIowACGlnFrAyMjJYM0amaskN4yVJP/3TDh2uo008QgthBEDIau90y9XaKYk1I5dq8TVjdFV+mk63d+n+N8rU3uk2uyTAizACIGR5DshzRNmUHBdtcjXWZo+y6V/vvFwp8dHadrRR//qXfWaXBHgRRgCErLqeKZq0xBhFRdlMrsb6cofG68nbp0qSVpQe1IYDdSZXBHQjjAAIWZ7Fq5xL4z8Lp+bo767Mk2FIP/rjVjX03GPATIQRACHrTI8RFq/60z/dMkmjhyWqxtWmH6/cxnZfmI4wAiBkcS5NYCTEOPTMXTMUbbdpza4a/f7zCrNLQoQjjAAIWXU9J/bSCt7/pgxP0Y9vmiBJeuLPu3SgtsnkihDJCCMAQhYjI4H13a/ka+64DLV2uLV0Fd1ZYR7CCICQxbk0gRUVZdM/3zFNCTF2bTrcQHdWmIYwAiBk1XNib8DlDo3XAwvGSZKWrd6tk6fZXYPgI4wACFlM0wTHd6/O17jMIapvbte/fLDX7HIQgQgjAEKSYRiqb2ZrbzBE26P0i9umSJJ+/3mFtlaeNLcgRBzCCICQdLq9S60d3eenMDISeFeNTtftM4bLMKR/fHu7utwsZkXwEEYAhCTPFE18tF0JMQ6Tq4kMSxdNVFKcQzuqXPr3z46YXQ4iiM9hZP369br11luVm5srm82mt99++6KvKS0tVWFhoeLi4jR69GitWLFiMLUCiCDeHiOMigTNsKRYPXxjgSTpXz7Yq9qmVpMrQqTwOYw0Nzdr+vTpevbZZwf0/PLyci1atEhz585VWVmZHnnkEd1///1auXKlz8UCiBzexaucSxNU37xqlKYOT1FTa6eWrd5jdjmIED6PfS5cuFALFy4c8PNXrFihkSNHavny5ZKkiRMnavPmzXrqqad0xx13+Pr2ACKEp8dIOotXg8oeZdMvbpui2577RG+VVekbM/NUNCbd7LIQ5gK+ZmTjxo0qLi7uc+3GG2/U5s2b1dHR0e9r2tra5HK5+jwARJYzPUYYGQm26XlD9a2rRkqS/ulPO9Te6Ta5IoS7gIcRp9OprKysPteysrLU2dmpurq6fl+zbNkypaSkeB95eXmBLhNAiKljZMRU/6t4gtITY3Sg9pRe+rjc7HIQ5oKym8Zms/X52HP+wdnXPZYuXarGxkbvo7KyMuA1AggtnjUjtII3R0pCtB5ZNFGS9Mxf96vqZIvJFSGcBTyMZGdny+l09rlWW1srh8Oh9PT+5yFjY2OVnJzc5wEgstSzm8Z0t18xXFdelqaWji49/l87zS4HYSzgYaSoqEglJSV9rq1Zs0YzZ85UdHR0oN8egEWd2U3DNI1ZbDabnrhtihxRNn2ws0Yf7qkxuySEKZ/DyKlTp7R161Zt3bpVUvfW3a1bt6qiokJS9xTL3Xff7X3+4sWLdeTIES1ZskS7d+/Wyy+/rJdeekkPPfSQf74DAGGpjnNpQkJBdpK+d3W+JGnpqu060cxBevA/n8PI5s2bNWPGDM2YMUOStGTJEs2YMUM//elPJUnV1dXeYCJJ+fn5Wr16tdatW6fLL79cTzzxhJ555hm29QI4L7fbUMNpz5oRRkbM9sD14zRmWKJqXG36X//xpXfdH+AvNsMCf6pcLpdSUlLU2NjI+hEgAjQ0t2vGE93Tu/t+sVAxDk6uMNuuYy7d9twnau9067FbJ+k7X8k3uyRYwEB/f/M3HEDI8SxeTYmPJoiEiEm5yXq0Z3fNstV7tKOq0eSKEE74Ww4g5LBeJDTdXTRKN0zKUnuXW/f/oUzNbZ1ml4QwQRgBEHI4lyY02Ww2/eqOacpOjtOhumb97B22+4aLptb+O6IHC2EEQMjx9hhhW2/ISU2M0fK7LleUTfqPLUf1p61VZpeES9DY0qH/vXq3ipZ9qMoTp02rgzACIOQwTRPaZo9O133XjZMkPfrWDh2pbza5Iviqo8ut3204rPn/slYvrD+kU22dpgZLwgiAkMOJvaHv/uvGatZlqTrV1qn7/1DGYXoWYRiGSnbV6Mbl6/XYOzvVcLpDYzOH6JX/Nkv3XjvWtLoIIwBCDufShD6HPUrL75qhlPhofXm0Uf9/yV6zS8JF7Khq1Df/7TP999c269DxZqUnxuiJ26bo/Qfm6toJmec9Ly4YHKa9MwCcB2tGrGH40Hj98x3TtPj1Lfo/pYf0lTEZmjd+mNll4SzOxlb9ywd7tarsqAxDinFE6XtX5+sH88coOS40jmUhjAAIOfWsGbGMm6Zk69uzR+r1Tyu05I9f6r0H5mpYEiHSbG63obLKk3p3W7V+//kRtXZ0T6P9zfRcPXxTgUakJphcYV+EEQAhp76ZaRor+cebJ2lTeYP21jTpe7/bpP/9tamaMjzF7LIiTmeXW5+Xn9B7O5z6YKdTtU1t3s/NHJWqR2+eqBkjU02s8PwIIwBCSnunW40t3T0PmKaxhrhou37zzRm6/bkN2na0Ubf85mN9bcZwPXRjgYYPjTe7vLDW1tmlTw7U6f0dTpXsqlHD6TP9QpJiHbpuYqa+enmuri0wd03IxRBGAIQUzwF59iibUuJDYz4bFzc+K0nvPTBXT63Zqz9tPaa3yqr07vZqffcr+fr/rg2dtQnhYmvlSb228bBKdtaoqVcn3NSEaBVPytZNU7I1Z2y6Yh12E6scOMIIgJBS17OtNy0xRlFRofsvOZwrLy1Bv75rhr53db5++e5ufVZ+QitKD+qPmyt1/3Vj9a3ZoxRtZxPnYHV0ufX+Dqde/qRcZRUnvdczk2J105Rs3TQ5W1fmp8lhwXtMGAEQUmgFb33TRgzVG/9jtv66u1bL3tutg8eb9bP/2qXfbTyiH980QTdOzgrpKYNQ09Dcrt9/XqH/u/GInK5WSVK03aZbp+Xqm1eN1BUjUy0f3AkjAEKKZ1tvBg3PLM1ms+n6SVmaXzBMb2yq1PK/7FN5XbMWv75FhaNSdd+1YzW/YBih5AL2Opv0yiflequsSm09TeUyhsToW1eN0rdmj1RmUpzJFfoPYQRASPGMjKQxMhIWHPYofXv2KN02Y7j+T+lB/dtHh7TlSIO+8+omTcxJ1g/mj9HNU3Nkt/i/7P1pR1Wjnnxvjz4+UOe9NmV4sr4zJ1+3TM+xzDoQXxBGAIQUzqUJT0NiHfqH4gJ9e/YovfjRIf37ZxXaXe3S/X8o09Nr9up/XjNGt18xPCx/0Q5UY0uHnl6zV//30yNyG1KUTbpxcra+e3W+Zo5KDetRJMIIgJByvIlpmnCWlRynR2+epHuvHavfbTiiVzaU63D9aS1dtV3L/7JP3796tL551UglxkbOryfDMPT21ir98t093gXct0zL0Y9vmqC8tNBqThYokfP/NgBLqOlZoJedHD7z4TjX0IQYPXD9OH1/br7+8HmFXvyoXE5Xq365ereeXXtAdxeN0q3TczUuc0hYjwjsq2nSP729Q5+Vn5AkjR6WqCe+OkVfGZthcmXBRRgBEFKqG1skSTkphJFIkBjr0PfnjtbfF43S22VVWlF6SOV1zfrNhwf0mw8PaFR6gm6YmKXrJ2Vp5qhUS25b7U9zW6ee+et+vfRxuTrdhuKio/TD67rDWSROVRFGAISUGlf3MHUWYSSixDrsunPWSH29ME/v73DqP7dU6pOD9TpSf1ovflyuFz8uV2pCtK6dkKniSVmaO26YJadyDMPQ+zucevzPu1Td2D0KeMOkLP30lkkRMyXTH+v9PwkgbDW1duhUTzdJpmkikz3Kppun5ejmaTlqbuvUR/uPa82uGn24p1YNpzu06osqrfqiSjGOKF09NkPfmJmn6ydmhvSISXunW5+V16tkV43+sqtGx3pCyIjUeP38byZrwcQskys0H2EEQMhw9vyQTo5zWPJfvfCvxFiHbpqSo5um5Kizy63NRxpUsqtGJbtqVHHitD7cU6sP99Rq+NB4fXv2KN01K0+pIbIlvLGlQ+v21qpkV41K9x7v07I9Ptqu78/N173XjlVcdORNyfSHv+0AQoanu2Q2UzQ4i8Mepdmj0zV7dLr+8eaJ2ldzSn/aWqU/fF6hqpMt+uf392j5X/bpq5fn6p45l2lybuBODe7scut0R5da2rt0ur1Lp9s7e/63S4frmlWyq0afHqpXp9vwviZjSKyun5ipGyZl6StjMwghZyGMAAgZnjn07BROesX52Ww2FWQn6eGbJuj+BeP0X18e0+82HtaOKpf+uPmo/rj5qGZdlqp75lymGydnX9J5OIfrmvXBTqfW7KrRweOndLq9S+093VAvZmzmEN0wKUvXT8zSjLyhlm/ZHkiEEQAho8YTRpLpMYKBiYu2629n5unrhSP0RUWDXt1wRO9tr9amww3adLhB2clxml8wTJNzkzUpN1kTc5KVEHP+X32GYWhXtUsf7KzRBzuc2lvTdN7n2qNsSoi2Kz7GroQYuxJiHEpLjNG88Rm6YVK28jMSA/EthyXCCICQUe1iZASDY7PZVDgqTYWj0lRz80T9+2cV+v1n3QfLvbGpstfzpPyMRE3KSdbk3BRNzk3WhJwkHa47rQ92OvXBTqeONrR4n2+Psmn26DTdODlbV+Wna0icwxtAYh1RYd0DJZgIIwBChmcBKz1GcCmykuO05IbxuvfaMSrde1zbjjZq57FG7TzmUm1Tmw4db9ah483687bqfl8fFx2leeOG6cbJ2VowMVNDE0JjUWw4I4wACBnORrqvwn9iHXYVT85W8eRs77XjTW3aeaxRu6pd2nnMpV3HXCqva1ZynEPXT8xS8eRsXTN+mOJjWGAaTIQRACGD3TQItGFJsZpfkKn5BZneay3tXYq220K6V0m4I4wACAmtHV060dx9Yi8jIwgmRkHMRwwEEBJqe9rAxzqiNDQh2uRqAAQTYQRASOh9QB47FIDIQhgBEBI860WymKIBIg5hBEBIYFsvELkIIwBCgpOGZ0DEIowACAlOWsEDEYswAiAkcEgeELkIIwBCQg0Nz4CIRRgBYLout6Hapu4+IyxgBSIPYQSA6epOtanLbcgeZVPGENaMAJGGMALAdJ71IplJsbJH0fAMiDSEEQCmc/Z0X2W9CBCZCCMATHdmWy9hBIhEhBEApqtmJw0Q0QgjAExXQyt4IKIRRgCYzrOAlUPygMhEGAFgOs+5NDl0XwUiEmEEgKkMw2ABKxDhCCMATHXydIfaOt2SpEwOyQMiEmEEgKk8UzTpiTGKi7abXA0AMxBGAJjKyeJVIOIRRgCYqpptvUDEI4wAMJVnmiaLMAJELMIIAFN5zqXJYZoGiFiEEQCmcrraJNEKHohkhBEApuLEXgCEEQCmYgErAMIIANM0t3WqqbVTkpRNK3ggYhFGAJjGs5NmSKxDQ2IdJlcDwCyEEQCmqfGcScMUDRDRCCMATFPNAXkARBgBYCLPNA0jI0BkI4wAMI2TnTQARBgBYKJqDskDIMIIABPVuBgZAUAYAWAiRkYASIMMI88995zy8/MVFxenwsJCffTRR+d97rp162Sz2c557NmzZ9BFA7C+9k636k51n0vDyAgQ2XwOI2+++aYefPBBPfrooyorK9PcuXO1cOFCVVRUXPB1e/fuVXV1tfcxbty4QRcNwPpqm7pHRWLsUUpLjDG5GgBm8jmMPP300/re976n73//+5o4caKWL1+uvLw8Pf/88xd8XWZmprKzs70Pu90+6KIBWJ9nJ01WSqxsNpvJ1QAwk09hpL29XVu2bFFxcXGf68XFxdqwYcMFXztjxgzl5ORowYIFWrt27QWf29bWJpfL1ecBILx4eozkJHMmDRDpfAojdXV16urqUlZWVp/rWVlZcjqd/b4mJydHL7zwglauXKlVq1apoKBACxYs0Pr168/7PsuWLVNKSor3kZeX50uZACzgzMgI60WASDeok6nOHlI1DOO8w6wFBQUqKCjwflxUVKTKyko99dRTmjdvXr+vWbp0qZYsWeL92OVyEUiAMFNNwzMAPXwaGcnIyJDdbj9nFKS2tvac0ZILmT17tvbv33/ez8fGxio5ObnPA0B48baCZ1svEPF8CiMxMTEqLCxUSUlJn+slJSWaM2fOgL9OWVmZcnJyfHlrAGHGyYm9AHr4PE2zZMkS/f3f/71mzpypoqIivfDCC6qoqNDixYsldU+xVFVV6bXXXpMkLV++XJdddpkmT56s9vZ2vf7661q5cqVWrlzp3+8EgKUQRgB4+BxG7rzzTtXX1+vxxx9XdXW1pkyZotWrV2vUqFGSpOrq6j49R9rb2/XQQw+pqqpK8fHxmjx5st59910tWrTIf98FAEtxuw1vK3imaQDYDMMwzC7iYlwul1JSUtTY2Mj6ESAM1Da16spf/lVRNmnvLxYq2s7JFEA4Gujvb34CAAi6msbuNvDDkmIJIgAIIwCCr7qxRRJTNAC6EUYABJ13vQiLVwGIMALABGcantEKHgBhBIAJvK3gmaYBIMIIABN4D8ljmgaACCMATMDICIDeCCMAgsowDEZGAPRBGAEQVK7WTp1u75LEbhoA3QgjAILKM0UzNCFacdF2k6sBEAoIIwCCysmZNADOQhgBEFROT/dVpmgA9CCMAAgqZ8+5NCxeBeBBGAEQVE6X51wauq8C6EYYARBUnlbw2SmxJlcCIFQQRgAEldMbRhgZAdCNMAIgqNhNA+BshBEAQdPa0aWTpzsksZsGwBmEEQBB45miSYixKznOYXI1AEIFYQRA0HgXrybHyWazmVwNgFBBGAEQNDWe9SJM0QDohTACIGh6j4wAgAdhBEDQMDICoD+EEQBBU91zLg2t4AH0RhgBEDSVJzyH5NHwDMAZhBEAQdHc1qm9NU2SpCnDk02uBkAoIYwACIovKhrU5TY0IjVeOYyMAOiFMAIgKD4vPyFJuvKyNJMrARBqCCMAgsITRmblE0YA9EUYARBwbZ1d2lp5UpJ0JWEEwFkIIwACbkdVo9o63coYEqPRGYlmlwMgxBBGAATcZz1TNDNHpXEmDYBzEEYABNwm1osAuADCCICA6nIb2nykQRI7aQD0jzACIKD2OpvU1NqpIbEOTcxJMrscACGIMAIgoDYd7p6iuWJUqhx2fuQAOBc/GQAE1JlmZ6kmVwIgVBFGAASMYRj6vGdkZBbrRQCcB2EEQMAcqT+t401tirFHaXreULPLARCiCCMAAsYzKjJtRIriou0mVwMgVBFGAASMd70I/UUAXABhBEDAeHbS0OwMwIUQRgAERK2rVUfqT8tmkwpHsZMGwPkRRgAEhGe9yMTsZCXHRZtcDYBQRhgBEBCbWC8CYIAIIwAC4jPCCIABIowA8LvGlg7trWmSRLMzABdHGAHgd1uOnJBhSPkZiRqWFGt2OQBCHGEEgN99Xt4gSZrFeTQABoAwAsDvPi+vl8QUDYCBIYwA8KvWji5tr2qUJF2Vn25yNQCsgDACwK/KKk6qo8tQVnKs8tLizS4HgAUQRgD4lbcF/GVpstlsJlcDwAoIIwD8yhNG6C8CYKAIIwD8prPLrS1HPDtpCCMABoYwAsBvdh5z6XR7l5LjHCrISjK7HAAWQRgB4De914tERbFeBMDAEEYA+M3nPefRzGK9CAAfEEYA+IXbbfQZGQGAgSKMAPCLg8dPqeF0h+KiozR1eIrZ5QCwEMIIAL/4vGdU5PK8oYpx8KMFwMDxEwOAX2wq9/QXoQU8AN8QRgBcshpXqz452H043pWsFwHgI4fZBQCwtg0H6nT/G2WqO9WujCGxumLUULNLAmAxhBEAg+J2G/rt2gP617/sk9uQJmQn6blvXaGEGH6sAPANPzUA+OxEc7t+9OZWle47Lkn6xswRevyrUxQXbTe5MgBWNKg1I88995zy8/MVFxenwsJCffTRRxd8fmlpqQoLCxUXF6fRo0drxYoVgyoWgPm+qGjQLc98pNJ9xxXriNKvvj5Nv/r6dIIIgEHzOYy8+eabevDBB/Xoo4+qrKxMc+fO1cKFC1VRUdHv88vLy7Vo0SLNnTtXZWVleuSRR3T//fdr5cqVl1w8gOAxDEMvf1yub6zYqGONrcrPSNTb935F35iZZ3ZpACzOZhiG4csLrrrqKl1xxRV6/vnnvdcmTpyo2267TcuWLTvn+T/+8Y/1zjvvaPfu3d5rixcv1pdffqmNGzcO6D1dLpdSUlLU2Nio5ORkX8oF4AdNrR16+D+36b0dTknSzVNz9OQdU5UUF21yZQBC2UB/f/u0ZqS9vV1btmzRT37ykz7Xi4uLtWHDhn5fs3HjRhUXF/e5duONN+qll15SR0eHoqPP/WHW1tamtra2Pt9MIKzcclQ7jjUG5GsDVmEYUlunW20dXWrp6FJrR5daO9xq7exSS3uX2jrdqj/VJldrp6LtNj26aKLumXOZbDYOwgPgHz6Fkbq6OnV1dSkrK6vP9aysLDmdzn5f43Q6+31+Z2en6urqlJOTc85rli1bpp///Oe+lDYopfuO650vjwX8fYBwMHxovJ795gzNGJlqdikAwsygdtOc/S8iwzAu+K+k/p7f33WPpUuXasmSJd6PXS6X8vL8Py99w6Qs5aXF+/3rAlYT67ArLjpKcdH2Mw/HmY/jo+0alzWERaoAAsKnMJKRkSG73X7OKEhtbe05ox8e2dnZ/T7f4XAoPb3/ttGxsbGKjY31pbRBuXV6rm6dnhvw9wEAAOfn026amJgYFRYWqqSkpM/1kpISzZkzp9/XFBUVnfP8NWvWaObMmf2uFwEAAJHF5629S5Ys0YsvvqiXX35Zu3fv1o9+9CNVVFRo8eLFkrqnWO6++27v8xcvXqwjR45oyZIl2r17t15++WW99NJLeuihh/z3XQAAAMvyec3InXfeqfr6ej3++OOqrq7WlClTtHr1ao0aNUqSVF1d3afnSH5+vlavXq0f/ehH+u1vf6vc3Fw988wzuuOOO/z3XQAAAMvyuc+IGegzAgCA9Qz09/eg2sEDAAD4C2EEAACYijACAABMRRgBAACmIowAAABTEUYAAICpCCMAAMBUhBEAAGAqwggAADCVz+3gzeBpEutyuUyuBAAADJTn9/bFmr1bIow0NTVJkvLy8kyuBAAA+KqpqUkpKSnn/bwlzqZxu906duyYkpKSZLPZ/PZ1XS6X8vLyVFlZyZk3QcD9Di7ud3Bxv4OL+x18g7nnhmGoqalJubm5ioo6/8oQS4yMREVFacSIEQH7+snJyfxhDiLud3Bxv4OL+x1c3O/g8/WeX2hExIMFrAAAwFSEEQAAYKqIDiOxsbF67LHHFBsba3YpEYH7HVzc7+DifgcX9zv4AnnPLbGAFQAAhK+IHhkBAADmI4wAAABTEUYAAICpCCMAAMBUER1GnnvuOeXn5ysuLk6FhYX66KOPzC4pLKxfv1633nqrcnNzZbPZ9Pbbb/f5vGEY+tnPfqbc3FzFx8dr/vz52rlzpznFhoFly5Zp1qxZSkpKUmZmpm677Tbt3bu3z3O45/7z/PPPa9q0ad7GT0VFRXrvvfe8n+deB86yZctks9n04IMPeq9xv/3rZz/7mWw2W59Hdna29/OBut8RG0befPNNPfjgg3r00UdVVlamuXPnauHChaqoqDC7NMtrbm7W9OnT9eyzz/b7+V/96ld6+umn9eyzz2rTpk3Kzs7WDTfc4D2DCL4pLS3Vvffeq08//VQlJSXq7OxUcXGxmpubvc/hnvvPiBEj9OSTT2rz5s3avHmzrrvuOn31q1/1/kDmXgfGpk2b9MILL2jatGl9rnO//W/y5Mmqrq72PrZv3+79XMDutxGhrrzySmPx4sV9rk2YMMH4yU9+YlJF4UmS8dZbb3k/drvdRnZ2tvHkk096r7W2thopKSnGihUrTKgw/NTW1hqSjNLSUsMwuOfBkJqaarz44ovc6wBpamoyxo0bZ5SUlBjXXHON8cADDxiGwZ/tQHjssceM6dOn9/u5QN7viBwZaW9v15YtW1RcXNznenFxsTZs2GBSVZGhvLxcTqezz72PjY3VNddcw733k8bGRklSWlqaJO55IHV1demNN95Qc3OzioqKuNcBcu+99+rmm2/W9ddf3+c69zsw9u/fr9zcXOXn5+uuu+7SoUOHJAX2flvioDx/q6urU1dXl7Kysvpcz8rKktPpNKmqyOC5v/3d+yNHjphRUlgxDENLlizR1VdfrSlTpkjingfC9u3bVVRUpNbWVg0ZMkRvvfWWJk2a5P2BzL32nzfeeENffPGFNm3adM7n+LPtf1dddZVee+01jR8/XjU1NfrFL36hOXPmaOfOnQG93xEZRjxsNlufjw3DOOcaAoN7Hxj33Xeftm3bpo8//vicz3HP/aegoEBbt27VyZMntXLlSt1zzz0qLS31fp577R+VlZV64IEHtGbNGsXFxZ33edxv/1m4cKH3v6dOnaqioiKNGTNGv/vd7zR79mxJgbnfETlNk5GRIbvdfs4oSG1t7TmJD/7lWZXNvfe/H/7wh3rnnXe0du1ajRgxwnude+5/MTExGjt2rGbOnKlly5Zp+vTp+vWvf8299rMtW7aotrZWhYWFcjgccjgcKi0t1TPPPCOHw+G9p9zvwElMTNTUqVO1f//+gP75jsgwEhMTo8LCQpWUlPS5XlJSojlz5phUVWTIz89XdnZ2n3vf3t6u0tJS7v0gGYah++67T6tWrdKHH36o/Pz8Pp/nngeeYRhqa2vjXvvZggULtH37dm3dutX7mDlzpr71rW9p69atGj16NPc7wNra2rR7927l5OQE9s/3JS1/tbA33njDiI6ONl566SVj165dxoMPPmgkJiYahw8fNrs0y2tqajLKysqMsrIyQ5Lx9NNPG2VlZcaRI0cMwzCMJ5980khJSTFWrVplbN++3fi7v/s7Iycnx3C5XCZXbk0/+MEPjJSUFGPdunVGdXW193H69Gnvc7jn/rN06VJj/fr1Rnl5ubFt2zbjkUceMaKioow1a9YYhsG9DrTeu2kMg/vtb//wD/9grFu3zjh06JDx6aefGrfccouRlJTk/d0YqPsdsWHEMAzjt7/9rTFq1CgjJibGuOKKK7xbIXFp1q5da0g653HPPfcYhtG9Peyxxx4zsrOzjdjYWGPevHnG9u3bzS3awvq715KMV155xfsc7rn/fPe73/X+3Bg2bJixYMECbxAxDO51oJ0dRrjf/nXnnXcaOTk5RnR0tJGbm2vcfvvtxs6dO72fD9T9thmGYVza2AoAAMDgReSaEQAAEDoIIwAAwFSEEQAAYCrCCAAAMBVhBAAAmIowAgAATEUYAQAApiKMAAAAUxFGAACAqQgjAADAVIQRAABgKsIIAAAw1f8D4PBX64eJJv4AAAAASUVORK5CYII=", "text/plain": [ "
" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "plot.plot(classical_rdf)\n", "plot.show()" ] }, { "cell_type": "markdown", "id": "44895d7d-30db-4451-8243-66581c066c25", "metadata": {}, "source": [ "## Ring Polymer Molecular Dynamics\n", "\n", "RPMD is based on the path integral formulation of quantum mechanics, which says that a quantum mechanical trajectory is equivalent to summing over all possible classical trajectories. As a faster approximation, we randomly sample from the set of possible trajectories by simulating many copies of the system at once. As the number of copies goes to infinity, the RPMD result converges to the exact quantum result.\n", "\n", "Running a simulation with RPMD is very similar to running an ordinary simulation, except that we must tell it how many copies to simulate. For parahydrogen we can get a reasonably well converged result with 12 copies. Other systems often require more, which can make them very expensive to simulate." ] }, { "cell_type": "code", "execution_count": 6, "id": "4cc52e77-6eb8-41b1-ac8b-259a6ce45e5f", "metadata": {}, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAiMAAAGdCAYAAADAAnMpAAAAOnRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjEwLjMsIGh0dHBzOi8vbWF0cGxvdGxpYi5vcmcvZiW1igAAAAlwSFlzAAAPYQAAD2EBqD+naQAAQZxJREFUeJzt3XtcVHX+P/DXmRkYrjPK/SqCeMFLipAGZt6SFovWLqtl3+y6m5tdjK3vRu2v2/Zd+rbVWnnJVu2ya2bfUrPWTLooXksIzAQVAQVlEIGcGUAGZub8/hgYI0AZnOHM5fV8PM6jOHMO857PQ+Hl53wugiiKIoiIiIgkIpO6ACIiIvJsDCNEREQkKYYRIiIikhTDCBEREUmKYYSIiIgkxTBCREREkmIYISIiIkkxjBAREZGkFFIX0Bdmsxk1NTUIDAyEIAhSl0NERER9IIoi9Ho9oqKiIJP13v/hEmGkpqYGsbGxUpdBRERE/VBdXY2YmJheX3eJMBIYGAjA8mFUKpXE1RAREVFf6HQ6xMbGWn+P98YlwkjnoxmVSsUwQkRE5GIuNcSCA1iJiIhIUgwjREREJCmGESIiIpIUwwgRERFJimGEiIiIJMUwQkRERJJiGCEiIiJJMYwQERGRpBhGiIiISFIMI0RERCQpm8NIfn4+srKyEBUVBUEQsHnz5kveYzAY8PTTTyMuLg5KpRLDhg3D2rVr+1MvERERuRmb96Zpbm7G+PHjcc899+CWW27p0z3z5s3DmTNnsGbNGiQmJqKurg5Go9HmYomIiMj92BxGMjMzkZmZ2efrt23bhp07d6KiogJBQUEAgKFDh9r6tkRETq2qoQUfF1YjOW4w0ocFQ6mQS10Skctw+K69W7ZsQWpqKl5++WX861//gr+/P2688Ub89a9/ha+vr6PfnojI4Y7U6vBfq79HfZMBABCoVGBWUhh+MzYC00aEwdebwYToYhweRioqKrB79274+Phg06ZNqK+vx4MPPojGxsZex40YDAYYDAbr1zqdztFlEhH1y8Hqc1i49ntoz7cjLtgPLW0mnNUbsLm4BpuLa+DjJcP0EZZgMjMpDCofL6lLJnI6Dg8jZrMZgiBg3bp1UKvVAIDXXnsNt956K5YvX95j70hubi6ef/55R5dGRHRZvqtowH3vFaDJYMSE2EF4755JCPRR4Ieqn7Htp1psO1yLUz+fx7bDlv/3kguYkhiCP04bhskJwVKXT+Q0HD61NzIyEtHR0dYgAgBJSUkQRRGnTp3q8Z6cnBxotVrrUV1d7egyiYhssvPYWdz1zvdoMhhxVUIQ/n3/ZKj9vCCTCUgdGoS/3DAau/57Bj5/+Go8NCMRiWEBaDeJ2HH0LG775368uv0ojCaz1B+DyCk4PIxMmTIFNTU1aGpqsp47duwYZDIZYmJierxHqVRCpVJ1OYiInMW2nzS4/70DaG03Y8bIULx7zyQEKLt3NAuCgLHRajx+3Uh8lT0NX2Vfg9+lxEAUgTe/OY55q/ahurFFgk9A5FxsDiNNTU0oLi5GcXExAKCyshLFxcWoqqoCYOnVWLhwofX6BQsWIDg4GPfccw9KSkqQn5+PJ554Avfeey8HsBKRy9lUdAqLPyhCu0nE9eMiserOVPh49W2AamJYIP7+u/F48/ZkBCoV+KHqHOa8sQv/+VHj4KqJnJvNYaSgoADJyclITk4GAGRnZyM5ORnPPPMMAECj0ViDCQAEBAQgLy8P586dQ2pqKu644w5kZWXhjTfesNNHICIaGP/efxLZHx2EySzi1pQYvHF7MrwVtncwZ42PwtZHpyJ5yCDoW41Y/MEPePKTH9HSxvWXyDMJoiiKUhdxKTqdDmq1Glqtlo9siEgSb+eX429bjwAA7kqLw7NZYyCTCZf1PdtNZrz+VRmW7zgOUQSGhfrjzdsnYnQUf86Re+jr72/uTUNEdAmrd1VYg8iD04fhuRsvP4gAgJdchsevG4l1909GuEqJ8rPNmLt8D97dUwkX+Hcikd0wjBARXYTBaMIbX5cBALJnj8B//2YUBOHyg8gvpQ8LwRePXoNrk8LQZjLjuc9K8MTHP9r1PYicGcMIEdFFfFNaB12rEREqHyyekeiw9wny98Y/F6bi+RvHQC4T8HHhKXxf2eiw9yNyJgwjREQXsbHoNABgbnI05HZ4NHMxgiDgrvShmH9lLADg718e4eMa8ggMI0REvWhsbsO3R+oAADdPjB6w931k5nAoFTIcOPEzdhw9O2DvSyQVhhEiol58drAGRrOIsdEqjAgPHLD3jVD74O70oQCAl788CrOZvSPk3hhGiIh60fmI5ubknleLdqRF04YhUKlAqUaH/xziomjk3hhGiIh6UH62CQerz0EuE3DjhKgBf//B/t74/TUJAIDX8o6hnfvYkBtjGCEi6sGmHyy9ItNGhCIkQClJDfdeHY9gf29U1jfj48KeNxYlcgcMI0REv2I2i9jU8YjmpuSBG7j6awFKBR7smE78+ldlaG03SVYLkSMxjBAR/cr3Jxpx+tx5BCoVmD06XNJa7pg8BFFqH9TqWvHv/SclrYXIURhGiIh+ZeMPlkcic8ZF9nlHXkfx8ZJjybUjAADLvz0OfWu7pPUQOQLDCBHRL7S2m7D1UC2AgV1b5GJunhiNhFB//NzSjjW7K6Uuh8juGEaIiH5he8kZNBmMiBnsiyuHBkldDgBAIZfhT7NHAgBW76pEY3ObxBUR2RfDCBHRL3Q+orkpOdouO/PaS+bYCIyNVqHJYMTKHcelLofIrhhGiIg61OlbsausHoC0s2h6IpMJeDzD0jvy3r6T0GjPS1wRkf0wjBARddhSXAOTWcSE2EFICA2Qupxupo0IxaT4ILQZzXjj6zKpyyGyG4YRIqIOnWuL3OIkA1d/TRAE/Pd1lt6RjwpOobK+WeKKiOyDYYSICMDRWj0O1+jgJRdwwxUDv/x7X6UODcLMUWEwmUW8lndM6nKI7IJhhIgIwMYiy8DVGSPDMNjfW+JqLq5z7MhnB2twgr0j5AYYRojI45nMIjZ37tA7ceB36LXV6CgVpo0IBQD8X2G1xNUQXT6GESLyePvKG3BGZ4Da1wszRoVKXU6fzL8yFgDwceEpGLmjL7k4hhEi8nida4tkjY+EUiHt8u99NSspDIP9vHBGZ7BORyZyVQwjROTRmg1GfPGTZfn3m5Kd/xFNJ6VCbq13wwE+qiHXxjBCRB7ty8O1ON9uwtBgP0wcMkjqcmwy70pLGPmq9AwamgwSV0PUfwwjROTRNv5wYeCqIDjP8u99MSpChfExahjNonWNFCJXxDBCRB6rocmAveWW8RZzJzjnQmeX8rtUy0DWDQeqIYqixNUQ9Q/DCBF5rG+PnoVZBEZHqjAk2E/qcvrlxglRUCpkKKtrQnH1OanLIeoXhhEi8lhfl54BAFybFCZxJf2n8vHCnHGRACxLxBO5IoYRIvJIBqMJ+cfOAgBmJYVLXM3lmdfxqOazgzVoaTNKXA2R7RhGiMgj7a9oRHObCWGBSoyLVktdzmWZHB+EIUF+aDIY8cWhWqnLIbKZzWEkPz8fWVlZiIqKgiAI2Lx5c5/v3bNnDxQKBSZMmGDr2xIR2VXnI5pZSWGQyVxrFs2vyWQC5qV2rDlSwDVHyPXYHEaam5sxfvx4LFu2zKb7tFotFi5ciFmzZtn6lkREdiWKIr4urQMAzBrl2o9oOt2SEgNBAL6vbEQlN88jF2NzGMnMzMSLL76Im2++2ab7HnjgASxYsABpaWm2viURkV0dqdXj9Lnz8PGSYUpiiNTl2EWk2hfXDLfsq/MxN88jFzMgY0beeecdlJeX49lnn+3T9QaDATqdrstBRGQvX5VYHtFcnRgCX2/X2IumL7h5Hrkqh4eRsrIyPPnkk1i3bh0UCkWf7snNzYVarbYesbGxDq6SiDzJV0c6HtG4+CyaX+PmeeSqHBpGTCYTFixYgOeffx4jRozo8305OTnQarXWo7qaXY5EZB91+lYc7FgcbNYo111fpCfcPI9cVd+6KvpJr9ejoKAARUVFeOihhwAAZrMZoihCoVBg+/btmDlzZrf7lEollEqlI0sjIg/1TcfA1fExaoSpfCSuxv7mXRmDtXsqrZvnBQfwZyk5P4eGEZVKhUOHDnU5t2LFCnzzzTf4+OOPER8f78i3JyLq5qtS93xE06lz87yDp7TYVHQa909NkLokokuyOYw0NTXh+PHj1q8rKytRXFyMoKAgDBkyBDk5OTh9+jTef/99yGQyjB07tsv9YWFh8PHx6XaeiMjRWttN2H28c9VV93pE80u/S43FwVNabDhQjfuujne53YjJ89g8ZqSgoADJyclITk4GAGRnZyM5ORnPPPMMAECj0aCqqsq+VRIR2cHe8nq0tpsRpfbB6EiV1OU4DDfPI1cjiC6w57ROp4NarYZWq4VK5b4/QIjIsXI2HsL676tw51Vx+Otc9+6dfWxDMTYVncbtk4Yg9+ZxUpdDHqqvv7+5Nw0ReQRRFPHNkQtLwLs7bp5HroRhhIg8wk+ndTijM8DPW46rEoKlLsfhfrl53sYfTktdDtFFMYwQkUfI69gY75rhofDxcp9VV3sjkwm4O30oAGBVfjlXZCWnxjBCRB7hl7v0eorbJsUiyN8b1Y3n8Z9DGqnLIeoVwwgRuT2N9jwO1+ggCMAMN1t19WL8vBW4d8pQAMCKb8thNjv9fAXyUAwjROT2vu5Y6GzikMEI8bAVSe9MG4oApQJHz+jxdceePETOhmGEiNzeVx74iKaT2tcLd6bFAQCWf3scLrCaA3kghhEicmstbUbsLW8AAFzrpkvAX8q9U+KhVMhQXH0O+yoapC6HqBuGESJya7vK6tFmNCM2yBfDwwKkLkcSoYFKzL/Ssu7Iim/LJa6GqDuGESJya1+VWB7RXJsU7tF7tPzhmgQoZAJ2H6/nEvHkdBhGiMhtmc0ivj1qGbTpqY9oOsUM9sNvJ0QDAFZ8e/wSVxMNLIYRInJbxafOob6pDYFKBa4cGiR1OZL74/QECAKwveQMys7opS6HyIphhIjcVucjmmkjQ+Gt4I+7xLBAXDc6AgCwcgfHjpDz4N9OInJbneuLePojml96cMYwAMCnB2tQ3dgicTVEFgwjROSWTp87j6Nn9JDLBEwfGSp1OU7jiphBmDo8BCaziFX57B0h58AwQkRu6WDHjJHRkSoM8vOWthgns3hGIgDgo4JTqNO3SlwNEcMIEbmpkhodAGBMlEriSpzP5PggTBwyCG1GM9bsrpS6HCKGESJyT4drtACA0Qwj3QiCYO0d+fe+k9C2tEtcEXk6hhEickslGkvPyOhIhpGezBwVhlERgWhuM+G9fSekLoc8HMMIEbmd+iYDzugMEARgFMNIjwRBwIMdvSPv7KlES5tR4orIkzGMEJHbKe3oFRka7I8ApULiapzX9eMiMTTYDz+3tOOD76qkLoc8GMMIEbmdwzV8RNMXcpmARdMs6478c1cFDEaTxBWRp2IYISK30zmThoNXL+2midGIUPngjM6ATwpPS10OeSiGESJyO9bBqwwjl6RUyPGHaxIAAG/tLIfRZJa4IvJEDCNE5FbOt5lQcbYJADCGj2n65LZJsQjy90ZVYws+/1EjdTnkgRhGiMitHKnVwSwCIQHeCA1USl2OS/DzVuC+q+MBAMu/PQ6zWZS4IvI0DCNE5FYuPKJRQxAEiatxHXemxSFQqUBZXRO2d+x2TDRQGEaIyK2UcCZNv6h8vLAwPQ4AsGLHcYgie0do4DCMEJFb4eDV/rt3Sjx8veT48ZQWu8rqpS6HPAjDCBG5DZNZxBGNHgB7RvojOECJ2ycNAQAs+/a4xNWQJ2EYISK3UVnfjPPtJvh6yREf4i91OS7p99fEw0su4PvKRhw40Sh1OeQhbA4j+fn5yMrKQlRUFARBwObNmy96/caNGzF79myEhoZCpVIhLS0NX375ZX/rJSLqVecjmlGRgZDLOHi1PyLVvrg1JQaAZWYN0UCwOYw0Nzdj/PjxWLZsWZ+uz8/Px+zZs7F161YUFhZixowZyMrKQlFRkc3FEhFdDAev2seiacMgE4AdR8/ip9NaqcshD2DzDlKZmZnIzMzs8/VLly7t8vXf/vY3fPrpp/jss8+QnJxs69sTEfXqcI3lF+eYKLXElbi2uGB/ZI2PwqfFNVix4zhW3JEidUnk5gZ8zIjZbIZer0dQUFCv1xgMBuh0ui4HEdHFiKLIPWns6MHpiQCAL36qxfE6vcTVkLsb8DDy6quvorm5GfPmzev1mtzcXKjVausRGxs7gBUSkSs6qzegobkNMgEYGR4odTkub2REIDJGh0MUgRU7yqUuh9zcgIaR9evX47nnnsOGDRsQFhbW63U5OTnQarXWo7q6egCrJCJXdLhj8GpCaAB8veUSV+MeFs+w9I58WlyD6sYWiashdzZgYWTDhg2477778NFHH+Haa6+96LVKpRIqlarLQUR0MZ2PaMbwEY3djI8dhKnDQ2Ayi1iVz94RcpwBCSPr16/H3XffjQ8++ADXX3/9QLwlEXkYzqRxjM7ekY8KTqFO1ypxNeSubA4jTU1NKC4uRnFxMQCgsrISxcXFqKqqAmB5xLJw4ULr9evXr8fChQvx6quv4qqrrkJtbS1qa2uh1XK6GBHZD5eBd4zJ8UFIiRuMNqMZb+dXSF0OuSmbw0hBQQGSk5Ot03Kzs7ORnJyMZ555BgCg0WiswQQAVq1aBaPRiMWLFyMyMtJ6PProo3b6CETk6ZoMRpxoaAbAnhF7EwQBD8+09I68v/8kx46QQ9i8zsj06dMvupvju+++2+XrHTt22PoWREQ2OaLRQRSBCJUPggOUUpfjdqaNCEX6sGDsLW/AK9uP4vXbuEYU2Rf3piEil8dHNI4lCAKempMEQbDMrDlYfU7qksjNMIwQkcvj4FXHGxutxk3J0QCA/9laetEeciJbMYwQkcvr7BnhtF7HejxjJJQKGb6vbEReyRmpyyE3wjBCRC6t3WTGkVrLcuV8TONYUYN8cf/UeADAS18cQbvJLHFF5C4YRojIpVWcbUab0YwApQKxg/2kLsftLZo2DMH+3qiob8b676sufQNRHzCMEJFLK9FY1ixKigyETCZIXI37C/TxwpLZIwAAS78qg661XeKKyB0wjBCRS7uwDLxa4ko8x21XxmJYqD8am9uwkpvokR0wjBCRSzvMmTQDzksuQ05mEgBgze5KnD53XuKKyNUxjBCRyxJFkWuMSGRWUhiuSghCm9GMV748KnU55OIYRojIZWm0rTjX0g6FTMDw8ACpy/EogiDg6TmjAQCbik7j0CnuN0b9xzBCRC6rc7xIYlgAlAq5xNV4nnExasydEAUA+J+tJVwIjfqNYYSIXJZ1vAgf0Ujm8etGwlshw/6KRnxdWid1OeSiGEaIyGV1Tuvl4FXpxAz2w71TLAuh/e2LUi6ERv3CMEJELuvCMvCc1iulB2cMw2A/L1ScbcaHB6qlLodcEMMIEbkk7fl2VDdappSyZ0RaKh8vLLnWshDaP/KOQXueC6GRbRhGiMgllXb0ikQP8oXaz0viamjB5CHWhdBe/6pM6nLIxTCMEJFLKuHgVafiJZfh2awxAID39p1A2Rm9xBWRK2EYISKXdGG8CMOIs7hmRChmjw6HySzi+c841Zf6jmGEiFxSCZeBd0r/7/rR8FbIsPt4PbaXnJG6HHIRDCNE5HLajGaU1VkeA/AxjXMZEuyH30+1TPV98T8laG03SVwRuQKGESJyOWV1erSbRKh9vRA9yFfqcuhXHpyeiAiVD6obz2P1rgqpyyEXwDBCRC6nVGPpFUmKDIQgCBJXQ7/mr1QgZ84oAMDyb8tRw1196RIYRojI5RyvawIAjAwPlLgS6s2N46Nw5dDBON9uQu4XR6Quh5wcwwgRuZzys5YwMiyMO/U6K0EQ8GzWGAgC8NnBGnxf2Sh1SeTEGEaIyOVYw0gow4gzGxutxu2ThgAAnt1yGCYzp/pSzxhGiMiltJvMqGpoAcAw4goezxgJlY8CpRod1n9fJXU55KQYRojIpZxsaIHRLMLfW45wlVLqcugSgvy98aeMkQCAV7cfxbmWNokrImfEMEJELuWX40U4k8Y13DF5CEaGB+Lnlnb8I++Y1OWQE2IYISKXwvEirkchl+HZrNEAgH/tP4kjtTqJKyJnwzBCRC6lvK4ZAJAQ4i9xJWSL9MQQzBkXAbMIPLflMPetoS4YRojIpXBar+t6ak4SlAoZ9lc0ct8a6sLmMJKfn4+srCxERUVBEARs3rz5kvfs3LkTKSkp8PHxQUJCAt56663+1EpEHk4URT6mcWExg/3w+6kJAICXvjiCdpNZ4orIWdgcRpqbmzF+/HgsW7asT9dXVlZizpw5mDp1KoqKivDUU0/hkUcewSeffGJzsUTk2c42GaBvNUImAHHBflKXQ/2waPowhAR4o7K+GR98x6m+ZKGw9YbMzExkZmb2+fq33noLQ4YMwdKlSwEASUlJKCgowCuvvIJbbrnF1rcnIg9WcdYyXiQ2yA8+XnKJq6H+CFAqsOTaEfjL5p+w9KtjuGliNFQ+XlKXRRJz+JiRffv2ISMjo8u56667DgUFBWhvb+/xHoPBAJ1O1+UgIuIjGvdw25WxSAwLwM8t7Vj+7XGpyyEn4PAwUltbi/Dw8C7nwsPDYTQaUV9f3+M9ubm5UKvV1iM2NtbRZRKRC+icSTMslDNpXJlCLsNTHbv6vrPnBKobWySuiKQ2ILNpfr0wUeeUrt4WLMrJyYFWq7Ue1dXVDq+RiJwfe0bcx4yRYUgfFow2oxmvbD8qdTkkMYeHkYiICNTW1nY5V1dXB4VCgeDg4B7vUSqVUKlUXQ4iIk7rdR+CIOCpOUkQBODT4hocrD4ndUkkIYeHkbS0NOTl5XU5t337dqSmpsLLi4OWiKhvzreZcPrceQDsGXEXY6PVuCk5GgDwP/8p5UJoHszmMNLU1ITi4mIUFxcDsEzdLS4uRlWVZYpWTk4OFi5caL1+0aJFOHnyJLKzs1FaWoq1a9dizZo1ePzxx+3zCYjII1TWN0MUgcF+Xgjy95a6HLKTxzNGQqmQ4fsTXAjNk9kcRgoKCpCcnIzk5GQAQHZ2NpKTk/HMM88AADQajTWYAEB8fDy2bt2KHTt2YMKECfjrX/+KN954g9N6icgmHC/inqIG+eL+qfEAuBCaJ7N5nZHp06dftCvt3Xff7XZu2rRp+OGHH2x9KyIiK4YR9/XH6YnYcKDauhDaXelDpS6JBhj3piEil1DeseDZsDBO63U3nQuhAcDSr45B19rzGlTkvhhGiMgllNexZ8Sd3XZlLIaF+nMhNA/FMEJETs9sFlFRbwkjCQwjbsmyEFoSAC6E5okYRojI6dVoz6O13QwvuYDYwb5Sl0MOMnNUGNISuBCaJ2IYISKn1zleZGiwPxRy/thyV4Ig4OnrLyyE9nHhKalLogHCv9VE5PQ4XsRzjI1W447JQwAAj//fQfz18xIYOd3X7TGMEJHTu7AMPGfSeIIXbhyLR2YmAgDW7K7E3e8cwLmWNomrIkdiGCEip8c1RjyLTCYgO2MkVtwxEb5ecuw+Xo8bl+3B0Vq91KWRgzCMEJHTq+hcY4RhxKPMGReJT/6YjpjBvqhqbMFNK/Zg20+1l76RXA7DCBE5NV1rO+r0BgBAQigf03ia0VEqbHnoaqQlBKOlzYRF/y7E61+VwWzmpnruhGGEiJxaZ69IuEqJQB/u9O2Jgvy98f59k3B3xzLx//jqGB5c9wOaDUZpCyO7YRghIqfGmTQEAF5yGZ67cQxevuUKeMtl2Ha4Fjev2MvF0dwEwwgROTUOXqVfmndlLNb/4SqEBipx9Iwe81ftQ2V9s9Rl0WViGCEip3YhjHC8CFmkxA3GloemYFioP2q0rZi/ah+O13GmjStjGCEip3Zht172jNAFkWpffPiHNIyKCESd3oD5q/ajVKOTuiyX1NpuwqYiaVe7ZRghIqfVbjLjZAOn9VLPQgOVWP/7qzA2WoWG5jbc/s/9OHRKK3VZLuVwjRZZb+7GYxsO4j8/aiSrg2GEiJxWdWML2k0i/LzliFD5SF0OOaHB/t5Yd/9VmBA7COda2rHgn/tRePJnqctyemaziLfzyzF3+R6U1TUhNFAJla9CsnoYRojIaXU+okkI9YdMJkhcDTkrta8X/n3/ZEwaGgS9wYiFa77DdxUNUpfltGrOnccdq7/D37YeQbtJxOzR4dj26FRMHR4qWU0MI0TktDoHryaE8BENXVyAUoF3770SVyeGoLnNhLve+R67y+qlLsvpfP5jDX6zNB/7Khrg6yVH7s3j8PadKQgOUEpaF8MIETktrjFCtvDzVmD1XamYMTIUre1m3PveAXxz5IzUZTkFfWs7sj8qxkMfFEHXasT4GDW2PjoVt08aAkGQvteRYYSInBZ36yVb+XjJ8dadKcgYHY42oxkP/KsQHxVUQxQ9d/n4ghONmPPGLmz84TRkAvDwzER8/Md0xIc4z98rhhEickqiKF6Y1sueEbKBUiHH8jsmImt8FNpNIv774x/xX2u+87jF0bQt7cjdWop5q/ahuvE8Ygb74qMH0vCnjJHwkjvXr3/phs4SEV1EQ3MbtOfbIQhwqn/BkWvwksuwdP4EjIoIxBtfl2HP8QZctzQfD89IxAPThsFb4Vy/jO1J39qOtbtPYPXuCuhbLfv33DwxGs/fOMZp93diGCEip9Q5XiRmsC98vOQSV0OuSC4TsHhGIm64IhJ/2fwTdpXV49W8Y/j0YA3+dtM4TIoPkrpEu2ppM+L9fSfx1s5ynGtpBwCMDA/EE9eNxLWjwyWu7uIYRojIKfERDdlLXLA/3r93ErYcrMFfPy/B8bomzFu1D/NTY5EzZxQG+XlLXeJlaW03Yf33VVj+bTnqmwwALNPhH7t2BK4fF+kS0+IZRojIKVVwgzyyI0EQ8NsJ0Zg+IgwvbTuC9d9XYUNBNb4qPYO/3JCEuROinWJWiS3ajGb8X2E1ln1zHBptKwAgNsgXS2aNwG8nREHhZONCLoZhhIicEnfrJUdQ+3kh9+ZxuGViNJ7adAjHzjThsQ0HsWpnBaaNCMXVw0Nw5dAgp3002G4y40BlI7aXnMG2n2pRq7OEkEi1Dx6eORy/S41xusGpfcEwQkRO6cJjGg5eJftLHRqEzx+eitW7K/D6V2U4UqvHkVo9VuVXQKmQ4cqhQbh6eAiuTgzB6EiVpI86WtqMyD92FtsPn8HXR+qgPd9ufS0kQImHZgzDbZOGOG2A6guGESJyOq3tJlT/3AKAu/WS43grZHhweiLmp8ZiV1k9dpXVY/fxszijM2D38XrsPm5ZwTXI3xvpw4IxfWQYMsaEQzUAM1Lqmwz45kgdth8+g11lZ2Ewmq2vBft749qkcGSMCceUxBCXDiGdGEaIyOmcaGiGKFr2HAn2d+3BheT8ggOUmJscjbnJ0R3r2zRZgklZPfZXNKCxuQ2f/6jB5z9q4L1Jhpkjw/DbCVGYMSrssoKA0WTGqZ/Po6K+CeV1zZb/nm1Gxdlm60DUTrFBvrhudAQyxkQgJW4w5C4wKNUWDCNE5HTK6y48onG1QYXk2gRBQGJYIBLDAnHPlHi0m8woqjqH3WVnsfWnWhyva8K2w7XYdrgWgUoFrhsbgd9OiEJaQnCvA0Z1re2oONuM8romlJ/tPJpxsqEZ7abeV4YdHanCdWMikDEmHKMiAt3670K/wsiKFSvw97//HRqNBmPGjMHSpUsxderUXq9ft24dXn75ZZSVlUGtVuM3v/kNXnnlFQQHB/e7cCJyXxy8Ss7CSy7DpPggTIoPwmOzR6BUo8enB0/js+Ia1Ghb8XHhKXxceAohAUrccEUk0oYFo+bceUvgqGtG+dkm1OkNvX5/pUKG+BB/DAsNwLBQfySEBiAh1B/xIf5Ou0CZI9gcRjZs2IAlS5ZgxYoVmDJlClatWoXMzEyUlJRgyJAh3a7fvXs3Fi5ciH/84x/IysrC6dOnsWjRItx///3YtGmTXT4EEbmXC3vSMIyQ8xAEAaOjVBgdpcKfrxuFgpM/49Pi0/jPIQ3qmwx4d+8JvLv3RI/3hgUqLYEjzB8JIQEYFmYJH1FqX5dYB8TRBNHG3YMmT56MiRMnYuXKldZzSUlJmDt3LnJzc7td/8orr2DlypUoLy+3nnvzzTfx8ssvo7q6uk/vqdPpoFarodVqoVKpbCmXiFzQDW/uwk+ndfjnwlTMdvKVI4najGbsPn4WnxbX4GitHkOC/JAYFtARPiw9HQMx6NUZ9fX3t009I21tbSgsLMSTTz7Z5XxGRgb27t3b4z3p6el4+umnsXXrVmRmZqKurg4ff/wxrr/++l7fx2AwwGC40K2l0+lsKZOIXJjZLFrHjCRwWi+5AG+FDDNHhWPmKAbn/rJpZZT6+nqYTCaEh3dt8PDwcNTW1vZ4T3p6OtatW4f58+fD29sbERERGDRoEN58881e3yc3Nxdqtdp6xMbG2lImEbmwWl0rzreboJAJGBLkJ3U5RDQA+rVM269H9Iqi2Oso35KSEjzyyCN45plnUFhYiG3btqGyshKLFi3q9fvn5ORAq9Vaj74+ziEi19c5XiQu2M8lV5IkItvZ9JgmJCQEcrm8Wy9IXV1dt96STrm5uZgyZQqeeOIJAMAVV1wBf39/TJ06FS+++CIiIyO73aNUKqFUKm0pjYjcROduvZxJQ+Q5bPpnh7e3N1JSUpCXl9flfF5eHtLT03u8p6WlBTJZ17eRyy2LxNg4dpaIPEDnMvAJDCNEHsPmPtDs7GysXr0aa9euRWlpKR577DFUVVVZH7vk5ORg4cKF1uuzsrKwceNGrFy5EhUVFdizZw8eeeQRTJo0CVFRUfb7JETkFo6d0QMARoQzjBB5CpvXGZk/fz4aGhrwwgsvQKPRYOzYsdi6dSvi4uIAABqNBlVVVdbr7777buj1eixbtgx/+tOfMGjQIMycORP/+7//a79PQURuQRTFX4SRQImrIaKBYvM6I1LgOiNEnqFO34pJ//M1ZAJQ8sJv3GIDMCJP1tff3xyqTkROo+xM50wafwYRIg/CMEJETuNoLceLEHkihhEichpldRwvQuSJGEaIyGlc6BlhGCHyJAwjROQURFG0jhkZGcEwQuRJGEaIyClotK3QG4xQyAQMDeYGeUSehGGEiJzC0Y71RRJC/eGt4I8mIk/Cv/FE5BTKOsLIcI4XIfI4DCNE5BSO1naMF2EYIfI4DCNE5BS4Jw2R52IYISLJmc0i1xgh8mAMI0QkueqfW9Daboa3QoY4zqQh8jgMI0QkuWMd64skhgZALhMkroaIBhrDCBFJrnO8CBc7I/JMDCNEJLnOZeCHc/AqkUdiGCEiyVl7Rjh4lcgjMYwQkaSMJjMqzjYD4EwaIk/FMEJEkjrR0II2kxl+3nJED/KVuhwikgDDCBFJ6tgvloGXcSYNkUdiGCEiSXUOXh0RxsGrRJ6KYYSIJNW58iqn9RJ5LoYRIpKUtWeEg1eJPBbDCBFJxmA04URDCwCGESJPxjBCRJKpONsMk1mEykeBcJVS6nKISCIMI0Qkmc6ZNCPCAyEInElD5KkYRohIMtYwwsGrRB6NYYSIJHO01rJbL6f1Enk2hhEikkzntF72jBB5NoYRIpJES5sRVY2WmTTcII/IszGMEJEkjtc1QRSBYH9vBAdwJg2RJ2MYISJJHDvTMV6EvSJEHq9fYWTFihWIj4+Hj48PUlJSsGvXrotebzAY8PTTTyMuLg5KpRLDhg3D2rVr+1UwEbmHC9N6OXiVyNMpbL1hw4YNWLJkCVasWIEpU6Zg1apVyMzMRElJCYYMGdLjPfPmzcOZM2ewZs0aJCYmoq6uDkaj8bKLJyLXxWm9RNTJ5jDy2muv4b777sP9998PAFi6dCm+/PJLrFy5Erm5ud2u37ZtG3bu3ImKigoEBQUBAIYOHXp5VRORyzvWsScNB68SkU2Padra2lBYWIiMjIwu5zMyMrB3794e79myZQtSU1Px8ssvIzo6GiNGjMDjjz+O8+fP9/o+BoMBOp2uy0FE7kPX2o4abSsAYDjDCJHHs6lnpL6+HiaTCeHh4V3Oh4eHo7a2tsd7KioqsHv3bvj4+GDTpk2or6/Hgw8+iMbGxl7HjeTm5uL555+3pTQiciFlHYNXI1Q+UPt6SVwNEUmtXwNYf72HhCiKve4rYTabIQgC1q1bh0mTJmHOnDl47bXX8O677/baO5KTkwOtVms9qqur+1MmETmpzvEiwzl4lYhgY89ISEgI5HJ5t16Qurq6br0lnSIjIxEdHQ21Wm09l5SUBFEUcerUKQwfPrzbPUqlEkol1x0gcledYYTjRYgIsLFnxNvbGykpKcjLy+tyPi8vD+np6T3eM2XKFNTU1KCpqcl67tixY5DJZIiJielHyUTk6jiThoh+yebHNNnZ2Vi9ejXWrl2L0tJSPPbYY6iqqsKiRYsAWB6xLFy40Hr9ggULEBwcjHvuuQclJSXIz8/HE088gXvvvRe+vr72+yRE5DKsG+SxZ4SI0I+pvfPnz0dDQwNeeOEFaDQajB07Flu3bkVcXBwAQKPRoKqqynp9QEAA8vLy8PDDDyM1NRXBwcGYN28eXnzxRft9CiJyGY3NbahvMgAAhnO3XiICIIiiKEpdxKXodDqo1WpotVqoVCqpyyGiy7C/ogG3vb0fMYN9sfvPM6Uuh4gcqK+/v7k3DRENqDIOXiWiX2EYIaIBddQ6rZdhhIgsGEaIaEAd6xi8OjKC40WIyIJhhIgGjCiKOFbXuVsve0aIyIJhhIgGzFm9Aeda2iETgGGh7BkhIguGESIaMMc69qQZGuwPHy+5xNUQkbNgGCGiAXOUe9IQUQ8YRohowByr5bReIuqOYYSIBgyn9RJRTxhGiGhAGIwmlGh0AIBx0epLXE1EnoRhhIgGRKlGjzajGYP9vBAX7Cd1OUTkRBhGiGhAFFf9DAAYHzsIgiBIXA0ROROGESIaEMXV5wAAybGDpS2EiJwOwwgRDYjOMDJhyCBJ6yAi58MwQkQO19jchhMNLQCACTGDpC2GiJwOwwgROdzBjl6RhBB/qP28pC2GiJwOwwgROVwRH9EQ0UUwjBCRw10YvDpI0jqIyDkxjBCRQ5nNonVa7wTOpCGiHjCMEJFDVTY0Q9dqhFIhw6hILgNPRN0xjBCRQxVXnQMAjI1Ww0vOHzlE1B1/MhCRQ3G8CBFdCsMIETkUFzsjokthGCEih2ltN6G0Y6feCewZIaJeMIwQkcP8dFoLo1lESIAS0YN8pS6HiJwUwwgROYx1vMgQ7tRLRL1jGCEih7GuvMpHNER0EQwjROQwndN6OZOGiC6GYYSIHKJO34rT585DEIBxMWqpyyEiJ8YwQkQO0dkrMjwsAIE+3KmXiHrHMEJEDnFhsTPuR0NEF9evMLJixQrEx8fDx8cHKSkp2LVrV5/u27NnDxQKBSZMmNCftyUiF1LU0TPCxc6I6FJsDiMbNmzAkiVL8PTTT6OoqAhTp05FZmYmqqqqLnqfVqvFwoULMWvWrH4XS0SuwWQW8eOpcwA4k4aILs3mMPLaa6/hvvvuw/3334+kpCQsXboUsbGxWLly5UXve+CBB7BgwQKkpaX1u1gicg3H65rQ3GaCn7ccI8K5Uy8RXZxNYaStrQ2FhYXIyMjocj4jIwN79+7t9b533nkH5eXlePbZZ/v0PgaDATqdrstBRK6juPpnAMAVMWrIZVzsjIguzqYwUl9fD5PJhPDw8C7nw8PDUVtb2+M9ZWVlePLJJ7Fu3TooFIo+vU9ubi7UarX1iI2NtaVMIpKYdbwIB68SUR/0awDrr5d1FkWxx6WeTSYTFixYgOeffx4jRozo8/fPycmBVqu1HtXV1f0pk4gkUsyVV4nIBn3rqugQEhICuVzerRekrq6uW28JAOj1ehQUFKCoqAgPPfQQAMBsNkMURSgUCmzfvh0zZ87sdp9SqYRSqbSlNCJyEs0GI46d0QOw7ElDRHQpNvWMeHt7IyUlBXl5eV3O5+XlIT09vdv1KpUKhw4dQnFxsfVYtGgRRo4cieLiYkyePPnyqicip/PjKS3MIhCl9kG4ykfqcojIBdjUMwIA2dnZuPPOO5Gamoq0tDS8/fbbqKqqwqJFiwBYHrGcPn0a77//PmQyGcaOHdvl/rCwMPj4+HQ7T0TuwfqIhr0iRNRHNoeR+fPno6GhAS+88AI0Gg3Gjh2LrVu3Ii4uDgCg0WguueYIEbmvoirLTBqOFyGivhJEURSlLuJSdDod1Go1tFotVCqV1OUQUS9EUcTkv32NOr0BHz2QhknxQVKXREQS6uvvb+5NQ0R2o9G2ok5vgFwmYFw0d+olor5hGCEiu+kcLzIqIhC+3nJpiyEil8EwQkR2w/EiRNQfDCNEZDdc7IyI+oNhhIjsot1kxqHTWgBc7IyIbMMwQkR2cbRWj9Z2MwJ9FEgICZC6HCJyIQwjRGQXRb94RCPjTr1EZAOGESKyi2LrTr2DJK2DiFwPwwgR2UVxNWfSEFH/MIwQ0WVraDKg/GwzAIYRIrIdwwgRXbY95Q0ALIudBQcoJa6GiFwNwwgRXbY9ZfUAgKsTQySuhIhcEcMIEV0WURSx+3hHGBnOMEJEtmMYIaLLcqKhBafPnYe3XMZdeomoXxhGiOiy7C47CwCYGDcIft4KiashIlfEMEJEl8X6iIbjRYionxhGiKjfjCYz9nbMpLl6eKjE1RCRq2IYIaJ+O3RaC32rESofBcZFq6Uuh4hcFMMIEfXb7o4pvenDQiDnfjRE1E8MI0TUb5zSS0T2wDBCRP3SbDDihyrLfjQcvEpEl4NhhIj65fvKRrSbRMQM9kVcsJ/U5RCRC2MYIaJ++eWUXkHgeBEi6j+GESLql87BqxwvQkSXi2GEiGxWp2/F0TN6CIJlJg0R0eVgGCEim+3peEQzJkqFIH9viashIlfHMEJENttdZll1dQpn0RCRHTCMEJFNRFHE7uOWzfGmJnIJeCK6fAwjRGST8rNNOKMzwFshQ+rQwVKXQ0RugGGEiGyyq2MWzaShQfDxkktcDRG5g36FkRUrViA+Ph4+Pj5ISUnBrl27er1248aNmD17NkJDQ6FSqZCWloYvv/yy3wUTkbQ6B69yvAgR2YvNYWTDhg1YsmQJnn76aRQVFWHq1KnIzMxEVVVVj9fn5+dj9uzZ2Lp1KwoLCzFjxgxkZWWhqKjososnooHVbjJjf0UjAGAq1xchIjsRRFEUbblh8uTJmDhxIlauXGk9l5SUhLlz5yI3N7dP32PMmDGYP38+nnnmmT5dr9PpoFarodVqoVKpbCmXiOyo4EQjbn1rHwb7eaHwL7Mh4069RHQRff39bVPPSFtbGwoLC5GRkdHlfEZGBvbu3dun72E2m6HX6xEUFGTLWxORE+gcL5KeGMIgQkR2o7Dl4vr6ephMJoSHh3c5Hx4ejtra2j59j1dffRXNzc2YN29er9cYDAYYDAbr1zqdzpYyichBOseLTOV4ESKyo34NYP31pliiKPZpo6z169fjueeew4YNGxAWFtbrdbm5uVCr1dYjNja2P2USkR3pW9tRVH0OAAevEpF92RRGQkJCIJfLu/WC1NXVdest+bUNGzbgvvvuw0cffYRrr732otfm5ORAq9Vaj+rqalvKJCIH+K6iESaziKHBfogN8pO6HCJyIzaFEW9vb6SkpCAvL6/L+by8PKSnp/d63/r163H33Xfjgw8+wPXXX3/J91EqlVCpVF0OIpLWbk7pJSIHsWnMCABkZ2fjzjvvRGpqKtLS0vD222+jqqoKixYtAmDp1Th9+jTef/99AJYgsnDhQrz++uu46qqrrL0qvr6+UKvVdvwoRORInWGEU3qJyN5sDiPz589HQ0MDXnjhBWg0GowdOxZbt25FXFwcAECj0XRZc2TVqlUwGo1YvHgxFi9ebD1/11134d133738T0BEDqfRnsfxuibIBCAtgWGEiOzL5nVGpMB1Roik9XHhKTz+fwcxPnYQPl08RepyiMhFOGSdESLyTLvLLLv0Xp0YLHElROSOGEaI6KJEUcTu4w0AgKsTQyWuhojcEcMIEV3Uj6e0qG8ywNdLjolxg6Quh4jcEMMIEV3Ue/tOAABmjw6HUiGXthgicksMI0TUqzp9Kz4/qAEA3Ht1vMTVEJG7Yhghol6t21+FNpMZyUMGYULsIKnLISI3xTBCRD0yGE1Y991JAMA9U9grQkSOwzBCRD36/KAG9U1tiFD5IHNshNTlEJEbYxghom5EUcTaPZUAgDvT4uAl548KInIc/oQhom4OnPgZh2t0UCpkWDBpiNTlEJGbYxghom7e6egVuXliNAb7e0tcDRG5O4YRIuqiurEFXx627K59dzoHrhKR4zGMEFEX/9p/EmYRuDoxBCMjAqUuh4g8AMMIEVk1G4z48PsqAMA9U4ZKWwwReQyGESKy2vjDKehajRga7IcZI8OkLoeIPATDCBEBAMxmEe/sPQEAuDt9KGQyQdqCiMhjMIwQEQAgv+wsKs42I1CpwK2psVKXQ0QehGGEiAAAa/ecAAD8LjUWAUqFtMUQkUdhGCEiHK9rQv6xsxAEyyMaIqKBxDBCRHh3r2WRs2uTwjEk2E/iaojI0zCMEHk4bUs7Pik8DYDTeYlIGgwjRB7uwwNVON9uwqiIQKQlBEtdDhF5IIYRIg9mNJnx/r6TAIB7p8RDEDidl4gGHsMIkYcymUX877YjOH3uPIL8vXHjhCipSyIiD8X5e0QeSNvSjkc3FGHH0bMAgEdnDYePl1ziqojIUzGMEHmYY2f0+MP7BTjR0AKlQoaXb70Cv50QLXVZROTBGEaIPMi2n2rxp4+K0dxmQvQgX6y6MwVjo9VSl0VEHo5hhMgDmM0iln5dhje+LgMAXJUQhOULJiI4QClxZUREDCNEbk/f2o7HNhTjq9I6AJa1RJ6akwQvOcevE5FzYBghcmPlZ5vwh/cLUH62Gd4KGf520zjcmhIjdVlERF0wjBC5IY32PDYVncbKb8uhNxgRqfbBW/+VgvGxg6QujYiom371065YsQLx8fHw8fFBSkoKdu3addHrd+7ciZSUFPj4+CAhIQFvvfVWv4olot61tBmxqegU/mv1d0h/6Ru8vO0o9AYjJg0NwpaHrmYQISKnZXPPyIYNG7BkyRKsWLECU6ZMwapVq5CZmYmSkhIMGTKk2/WVlZWYM2cOfv/73+Pf//439uzZgwcffBChoaG45ZZb7PIhiDyV2Szi+xON+KTwFLYe0qC5zWR9bVJ8EG6dGIObJkZzfAgROTVBFEXRlhsmT56MiRMnYuXKldZzSUlJmDt3LnJzc7td/+c//xlbtmxBaWmp9dyiRYtw8OBB7Nu3r0/vqdPpoFarodVqoVKpbCmXyC20m8zQnm+/cLS0o6j6HDb+cAqnfj5vvW5IkB9unhiNm5NjuPsuEUmur7+/beoZaWtrQ2FhIZ588sku5zMyMrB3794e79m3bx8yMjK6nLvuuuuwZs0atLe3w8vLq9s9BoMBBoOhy4dxhE8KT+GnGq1DvjdRT0QREEURRrMIk/nX/zWj3WT5uqXNCO15I7QtbdCeb+/S4/FrAUoFrh8XiVtSYnDl0MHcX4aIXI5NYaS+vh4mkwnh4eFdzoeHh6O2trbHe2pra3u83mg0or6+HpGRkd3uyc3NxfPPP29Laf2y89hZbDlY4/D3IbKXQB8FBvl5Qe3rhUi1L264IhIZoyPg682l3InIdfVrNs2v/+UliuJF/zXW0/U9ne+Uk5OD7Oxs69c6nQ6xsbH9KfWiZo8OR2yQr92/L9HFyAUBcpkMCrkAuUyAQvbL/8qgkAlQeskwyM8bal8vDPK1hA+VrxfkMvZ6EJH7sSmMhISEQC6Xd+sFqaur69b70SkiIqLH6xUKBYKDg3u8R6lUQql0/MqQWeOjkDWeO5USERFJyaYh9t7e3khJSUFeXl6X83l5eUhPT+/xnrS0tG7Xb9++HampqT2OFyEiIiLPYvN8v+zsbKxevRpr165FaWkpHnvsMVRVVWHRokUALI9YFi5caL1+0aJFOHnyJLKzs1FaWoq1a9dizZo1ePzxx+33KYiIiMhl2TxmZP78+WhoaMALL7wAjUaDsWPHYuvWrYiLiwMAaDQaVFVVWa+Pj4/H1q1b8dhjj2H58uWIiorCG2+8wTVGiIiICEA/1hmRAtcZISIicj19/f3NZRmJiIhIUgwjREREJCmGESIiIpIUwwgRERFJimGEiIiIJMUwQkRERJJiGCEiIiJJMYwQERGRpBhGiIiISFI2Lwcvhc5FYnU6ncSVEBERUV91/t6+1GLvLhFG9Ho9ACA2NlbiSoiIiMhWer0earW619ddYm8as9mMmpoaBAYGQhAEu31fnU6H2NhYVFdXc8+bAcD2Hlhs74HF9h5YbO+B1582F0URer0eUVFRkMl6HxniEj0jMpkMMTExDvv+KpWKf5gHENt7YLG9Bxbbe2CxvQeerW1+sR6RThzASkRERJJiGCEiIiJJeXQYUSqVePbZZ6FUKqUuxSOwvQcW23tgsb0HFtt74DmyzV1iACsRERG5L4/uGSEiIiLpMYwQERGRpBhGiIiISFIMI0RERCQpjw4jK1asQHx8PHx8fJCSkoJdu3ZJXZJbyM/PR1ZWFqKioiAIAjZv3tzldVEU8dxzzyEqKgq+vr6YPn06Dh8+LE2xbiA3NxdXXnklAgMDERYWhrlz5+Lo0aNdrmGb28/KlStxxRVXWBd+SktLwxdffGF9nW3tOLm5uRAEAUuWLLGeY3vb13PPPQdBELocERER1tcd1d4eG0Y2bNiAJUuW4Omnn0ZRURGmTp2KzMxMVFVVSV2ay2tubsb48eOxbNmyHl9/+eWX8dprr2HZsmU4cOAAIiIiMHv2bOseRGSbnTt3YvHixdi/fz/y8vJgNBqRkZGB5uZm6zVsc/uJiYnBSy+9hIKCAhQUFGDmzJn47W9/a/2BzLZ2jAMHDuDtt9/GFVdc0eU829v+xowZA41GYz0OHTpkfc1h7S16qEmTJomLFi3qcm7UqFHik08+KVFF7gmAuGnTJuvXZrNZjIiIEF966SXrudbWVlGtVotvvfWWBBW6n7q6OhGAuHPnTlEU2eYDYfDgweLq1avZ1g6i1+vF4cOHi3l5eeK0adPERx99VBRF/tl2hGeffVYcP358j685sr09smekra0NhYWFyMjI6HI+IyMDe/fulagqz1BZWYna2touba9UKjFt2jS2vZ1otVoAQFBQEAC2uSOZTCZ8+OGHaG5uRlpaGtvaQRYvXozrr78e1157bZfzbG/HKCsrQ1RUFOLj43HbbbehoqICgGPb2yU2yrO3+vp6mEwmhIeHdzkfHh6O2tpaiaryDJ3t21Pbnzx5UoqS3IooisjOzsbVV1+NsWPHAmCbO8KhQ4eQlpaG1tZWBAQEYNOmTRg9erT1BzLb2n4+/PBD/PDDDzhw4EC31/hn2/4mT56M999/HyNGjMCZM2fw4osvIj09HYcPH3Zoe3tkGOkkCEKXr0VR7HaOHINt7xgPPfQQfvzxR+zevbvba2xz+xk5ciSKi4tx7tw5fPLJJ7jrrruwc+dO6+tsa/uorq7Go48+iu3bt8PHx6fX69je9pOZmWn9/3HjxiEtLQ3Dhg3De++9h6uuugqAY9rbIx/ThISEQC6Xd+sFqaur65b4yL46R2Wz7e3v4YcfxpYtW/Dtt98iJibGep5tbn/e3t5ITExEamoqcnNzMX78eLz++utsazsrLCxEXV0dUlJSoFAooFAosHPnTrzxxhtQKBTWNmV7O46/vz/GjRuHsrIyh/759sgw4u3tjZSUFOTl5XU5n5eXh/T0dImq8gzx8fGIiIjo0vZtbW3YuXMn276fRFHEQw89hI0bN+Kbb75BfHx8l9fZ5o4niiIMBgPb2s5mzZqFQ4cOobi42HqkpqbijjvuQHFxMRISEtjeDmYwGFBaWorIyEjH/vm+rOGvLuzDDz8Uvby8xDVr1oglJSXikiVLRH9/f/HEiRNSl+by9Hq9WFRUJBYVFYkAxNdee00sKioST548KYqiKL700kuiWq0WN27cKB46dEi8/fbbxcjISFGn00lcuWv64x//KKrVanHHjh2iRqOxHi0tLdZr2Ob2k5OTI+bn54uVlZXijz/+KD711FOiTCYTt2/fLooi29rRfjmbRhTZ3vb2pz/9SdyxY4dYUVEh7t+/X7zhhhvEwMBA6+9GR7W3x4YRURTF5cuXi3FxcaK3t7c4ceJE61RIujzffvutCKDbcdddd4miaJke9uyzz4oRERGiUqkUr7nmGvHQoUPSFu3CemprAOI777xjvYZtbj/33nuv9edGaGioOGvWLGsQEUW2taP9Ooywve1r/vz5YmRkpOjl5SVGRUWJN998s3j48GHr645qb0EURfHy+laIiIiI+s8jx4wQERGR82AYISIiIkkxjBAREZGkGEaIiIhIUgwjREREJCmGESIiIpIUwwgRERFJimGEiIiIJMUwQkRERJJiGCEiIiJJMYwQERGRpBhGiIiISFL/HzgPu3GjE9RqAAAAAElFTkSuQmCC", "text/plain": [ "
" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "copies = 12\n", "integrator = RPMDIntegrator(copies, temperature, 1.0/picosecond, 1*femtosecond)\n", "context = Context(system, integrator)\n", "context.setPositions(positions)\n", "LocalEnergyMinimizer.minimize(context)\n", "context.setVelocitiesToTemperature(temperature)\n", "integrator.step(1000) # Equilibrate before collecting data\n", "rpmd_rdf = compute_rdf(context)\n", "plot.plot(rpmd_rdf)\n", "plot.show()" ] }, { "cell_type": "markdown", "id": "2faa0c92-139c-48fa-9def-bb1bc2c7396d", "metadata": {}, "source": [ "We can see at a glance that the distribution is significantly different from the classical one with a lower, broader peak. The simulation takes much longer to run, since it must repeat all the calculations for every copy in each time step." ] }, { "cell_type": "markdown", "id": "96ebb0f5-fdb0-4c8e-98ee-45584e94e072", "metadata": {}, "source": [ "## Adaptive Quantum Thermal Bath\n", "\n", "adQTB is a much faster alternative based on a different approximation (https://doi.org/10.1021/acs.jctc.8b01164). It uses a colored (i.e. frequency dependent) noise thermostat whose energy at each frequency matches the average energy of a quantum harmonic oscillator:\n", "\n", "$\\theta(\\omega, T) = \\hbar \\omega \\left( \\frac{1}{2} + \\frac{1}{e^{\\hbar \\omega / k_B T} - 1} \\right)$\n", "\n", "For harmonic degrees of freedom this produces the correct energy distribution. For non-harmonic ones it is an approximation, but still often a very good one.\n", "\n", "A problem that can happen with this approach is zero point energy leakage. The thermostat drives the system toward the desired quantum energy distribution, but the classical dynamics of the system causes it to continuously relax toward the\n", "classical distribution. The result is a distribution in between the quantum and classical ones with too much energy in the low frequency modes and too little energy in the high frequency modes.\n", "\n", "The adaptive version of the algorithm compensates for this by monitoring the energy distribution in the system and continuously adjusting the amount of noise at each frequency to bring the system closer to the intended distribution. In practice this works as follows.\n", "\n", "1. Divide the trajectory into short segments, typically on the order of 1000 time steps.\n", "2. During each segment, monitor the distribution of energy across modes.\n", "3. At the end of each segment, adjust the amount of noise at each frequency accordingly.\n", "\n", "This procedure tends to be very effective, but it does make adQTB slightly more complicated to use. Let's start by creating the integrator." ] }, { "cell_type": "code", "execution_count": 7, "id": "37a5a1ae-d772-4d04-b026-69bd6e14766a", "metadata": {}, "outputs": [], "source": [ "integrator = QTBIntegrator(temperature, 20.0/picosecond, 1*femtosecond)" ] }, { "cell_type": "markdown", "id": "5b1aec90-f68f-4f0e-880c-53f93a8f720e", "metadata": {}, "source": [ "Notice that we are using a larger friction coefficient than we did for the other integrators. This is usually required for adQTB. If the friction is too low, it may be impossible to fully compensate for zero point energy leakage. Even reducing the noise magnitude to zero would leave too much energy in some modes. This can slow down motions and cause sampling to take longer, but it still ends up being much faster than RPMD.\n", "\n", "Next we tell it how long the individual segments of the trajectory should be. In this case we set them to 0.5 ps." ] }, { "cell_type": "code", "execution_count": 8, "id": "b0684b95-6187-4244-8026-8bad564d8b8c", "metadata": {}, "outputs": [], "source": [ "integrator.setSegmentLength(0.5*picosecond)" ] }, { "cell_type": "markdown", "id": "bc7e7e9d-5be9-4324-977f-9985f8a07c88", "metadata": {}, "source": [ "By default, the adaptation algorithm keeps track of a different noise spectrum for every particle in the system and adapts each one separately. When multiple particles are expected to behave identically, you should tell it to use the same spectrum for all of them. This allows it to average the data from them, making the algorithm more stable and allowing adaptation to happen more quickly.\n", "\n", "To do this, we call `setParticleType()` to give every particle a type. If two particles have the same type, they share the same noise spectrum. Since every particle in our system is identical, we give them all the same type." ] }, { "cell_type": "code", "execution_count": 9, "id": "ccd82573-f5b2-46b0-8102-6a8f1e623027", "metadata": {}, "outputs": [], "source": [ "for i in range(particles):\n", " integrator.setParticleType(i, 0)" ] }, { "cell_type": "markdown", "id": "b2edaad1-48ae-4f8b-870b-3c4c11ed2ef7", "metadata": {}, "source": [ "Finally we must tell it how quickly to adjust the noise spectrum. We do this by setting the adaptation rate. Larger values lead to faster adaptation, but too large a value can lead to noisy results or even prevent it from converging. There is no easy rule for picking the adaptation rate. Usually you need to select it through trial and error." ] }, { "cell_type": "code", "execution_count": 10, "id": "e72d5c06-44b9-40e5-b064-0bd8d66d4cb0", "metadata": {}, "outputs": [], "source": [ "integrator.setDefaultAdaptationRate(0.5)" ] }, { "cell_type": "markdown", "id": "6ef96491-98b2-4a2c-9a47-7a422b078ea3", "metadata": {}, "source": [ "It also is possible to set different adaptation rates for different particle types by calling `setTypeAdaptationRate()`. For example, if there are many more particles of one type than other, you might want a larger adaptation rate for the type with more particles. In most cases, though, it is fine to use the same adaptation rate for all types. Since our system only has a single particle type, this is not relevant to it anyway.\n", "\n", "Now that our integrator is set up, let's create the Context. This part is exactly the same as for the other integrators." ] }, { "cell_type": "code", "execution_count": 11, "id": "e554caa7-f01e-40b7-87c3-8843ebae19be", "metadata": {}, "outputs": [], "source": [ "context = Context(system, integrator)\n", "context.setPositions(positions)\n", "LocalEnergyMinimizer.minimize(context)\n", "context.setVelocitiesToTemperature(temperature)" ] }, { "cell_type": "markdown", "id": "b1e8404e-9e02-4391-b990-280032997b29", "metadata": {}, "source": [ "Now we need to equilibrate the system. In the earlier cases we just ran a short simulation to let the particle positions and velocities relax. For adQTB, equilibration serves another important function: it gives the adaptation algorithm time to select the noise spectra. This may require a much longer equilibration, and you should monitor the progress of it to make sure the result is truly converged before you start collecting data." ] }, { "cell_type": "code", "execution_count": 12, "id": "6929bc65-969c-4874-8f85-a769389cd0fd", "metadata": {}, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAjoAAAGfCAYAAABWcXgAAAAAOnRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjEwLjMsIGh0dHBzOi8vbWF0cGxvdGxpYi5vcmcvZiW1igAAAAlwSFlzAAAPYQAAD2EBqD+naQAAjwxJREFUeJzs3XecXHW9//HXOdNnd2a29+ym90p6CE0gEBHB7hUp96qoIIL8rgW93itXrxHvvaiIol6VXPVSVKoNCEISAgHSSW+7m02215mdXs75/fGdnc0mm5ANIZvMfp6PRx7Jzs6c+Z7ZyZ73fL5NM03TRAghhBAiC+nD3QAhhBBCiHeLBB0hhBBCZC0JOkIIIYTIWhJ0hBBCCJG1JOgIIYQQImtJ0BFCCCFE1pKgI4QQQoisJUFHCCGEEFlLgo4QQgghspYEHSGEEEJkLetQ7rxixQqefPJJ9uzZg8vlYsmSJdx3331MmjQpcx/TNLn33nv5xS9+QXd3NwsXLuQnP/kJ06ZNO+FxV65cyT/+4z8ed3skEsHpdJ5S2wzDoKmpCY/Hg6ZpQzktIYQQQgwT0zTp7e2loqICXT/z9ZchBZ01a9Zw++23M3/+fJLJJN/4xjdYtmwZu3btIicnB4Dvf//73H///axcuZKJEyfyne98hyuvvJK9e/fi8XhOeGyv18vevXsH3HaqIQegqamJUaNGDeV0hBBCCHGOOHz4MFVVVWf8uNo72dSzvb2dkpIS1qxZw8UXX4xpmlRUVHDXXXfx1a9+FYBYLEZpaSn33Xcfn/3sZwc9zsqVK7nrrrvo6ek53abg9/vJy8vj8OHDeL3e0z6OEEIIIc6eQCDAqFGj6OnpwefznfHjD6micyy/3w9AQUEBAHV1dbS0tLBs2bLMfRwOB5dccgmvvfbaCYMOQDAYpKamhlQqxezZs/n2t7/NnDlzTnj/WCxGLBbLfN3b2wuoypAEHSGEEOL88m4NOzntzjDTNLn77rtZunQp06dPB6ClpQWA0tLSAfctLS3NfG8wkydPZuXKlTz77LM8+uijOJ1OLrzwQvbv33/Cx6xYsQKfz5f5I91WQgghhDjWaQedL3zhC7z11ls8+uijx33v2FRmmuZJk9qiRYv45Cc/yaxZs7jooov4/e9/z8SJE/nxj398wsfcc889+P3+zJ/Dhw+f7qkIIYQQIkudVtfVHXfcwbPPPsvatWsHDBwqKysDVGWnvLw8c3tbW9txVZ6T0XWd+fPnn7Si43A4cDgcp9F6IYQQQowUQ6romKbJF77wBZ588kleeuklxowZM+D7Y8aMoaysjFWrVmVui8fjrFmzhiVLlgzpebZu3TogLAkhhBBCDNWQKjq33347jzzyCM888wwejycz7sbn8+FyudA0jbvuuovvfve7TJgwgQkTJvDd734Xt9vNJz7xicxxbrrpJiorK1mxYgUA9957L4sWLWLChAkEAgEeeOABtm7dyk9+8pMzeKpCCCGEGGmGFHQeeughAC699NIBtz/88MPccsstAHzlK18hEolw2223ZRYMfOGFFwasodPQ0DBgUaCenh5uvfVWWlpa8Pl8zJkzh7Vr17JgwYLTPC0hhBBCiHe4js65JBAI4PP58Pv9Mr1cCCGEOE+829dv2etKCCGEEFlLgo4QQgghspYEHSGEEEJkLQk6QgghhMhaEnSEEEIIkbUk6LyN3tWrCTz/wnA3QwghhBCn4R3tXp7tonv2cOTzt4Guk/PqOix5ecPdJCGEEEIMgVR0TsA0TVrvuw9ME1Ipkh0dw90kIYQQQgyRBJ0TCL68mvD61zNfp7q7h7E1QgghhDgdEnQGYcbjtH3/+wNuS0rQEUIIIc47EnQG4f/LX4nX12MpLMSd3m8r1d0zvI0SQgghxJDJYORBxOvrAfBedRXBsEnEeVC6roQQQojzkASdQfR2tQDQbjV4o3cRqbkXUNq5a5hbJYQQQoihkq6rQdQ3q1CzoUsjlrKRtLnxd8eHuVVCCCGEGCoJOoOI+bsx0QhHpmZu6w2khrFFQgghhDgdEnQGkQr20lUwBS1VkLktFLEMY4uEEEIIcTok6BwjkoxgCcc4UnkJAHr6FQol7MPYKiGEEEKcDgk6x6j112Ixi+ksnI6JwZQLPACETdcwt0wIIYQQQyVB5xgHew4S9lwAQJt7N94JqpITseVhxk88IHn/xlbq3pJtIoQQQohziQSdYxzo3k/CUQ5Aa+5BYoUJAKLOAhInWEvH3x7mhV/u5Pn/2UEqZZy1tgohhBDi5CToHONQ+36irhIAOnJbCTr9AKSsLiItA4NOrb+Wr679Kts2H1T3SRhEexNnt8FCCCGEOCEJOsdobj5A2FUKQLunjc5kO/ZUCAD/ka4B9/3Nzt/w17q/8vqG7ZnbwgFZb0cIIYQ4V0jQOUo4ESbaHiFldYKZwu/spC3chpt00GkJDrh/nb8O3dCxNfv6j9ErQUcIIYQ4V0jQOUqtv5bCkOq2csa7MfQUHZEOcqwxAHq7IgPuXx+opyRYgz3VPyMrIhUdIYQQ4pwhQeco+7v3kx9VQceVUt1U7ZF2cpxqgHFvTzJzX3/MT1e0iyr/5AHHkIqOEEIIce6QoHOUgz0H8cRU0HETAKAj0oHHo16mYKj/vvWBegCq/VMAiFh71d9S0RFCCCHOGRJ0jnLAfwB3Qg1EzrWp7qq2cBu5+WotnVC8f7P3en899qSL4t5R6rFFmwGp6AghhBDnEgk6RznYcxCnoSo6vhy1iWd3tJvcIicA4ZQT0zQBVdGp9E9EQ8fmDHLB3gZAKjpCCCHEuUSCTlrKSNHR24nFVBt5FnptWDQLJiYUaOo+mo1oSK2TU+evY8aRMQCUHtjI1EOqqyvQHRrk6EIIIYQYDhJ00gLxAJ5oEZqmY0lGyPE5KHQVAhDMNbHHetT9OqKA6rqq6CkGwGsJohlqjE6wJ3z2Gy+EEEKIQUnQSeuJ9ZAXSQ9EDrdi8Xgodqkg0+VI4Ip2AtDbGSVlpGjobcBqqiDUXp1Lq0t1dSWjOqZhDsMZCCGEEOJYEnTSBgadNnRPbibotNkiOKNqurm/tZemYBOJVJKUTQUdu7mLVrsbAA09070lhBBCiOElQSetO9pNXnoNHXckXdFxq6DTqvfijKmg09MYoC5QR1FvLqZuRzNTTMg7xOipM7Am1MrJ2+pqh+ckhBBCCDGABJ20oys6OeFW9Nz+ik57pAOPpsbgdDeHqPPXUdNRBIAz1sVEZxOLllyAPa7u8/ib64fhDIQQQghxLAk6ad3R7gFjdPQ9f6DIVC9PR6QDn11t/9DdEae+p45yv+q2ykm2Y0sGsRfnZoJOXWs9/rB0XwkhhBDDTYJOWk8wgCOlxtk4o11Y9v6eku1PAGobCG8uYBrE49B4+CAFIRV0nKhByvacOPaEmmLuNnv5w6bDZ/8khBBCCDGABJ20QLeaFq6nYlhSUXSbSdGRLQC0h9ux53txR9oB8LclyI2roGO1qXVzrHRiT6p/F0YcrGvYcbZPQQghhBDHkKCTFgmoriZHrAcN0G0GxSk1Zbwz2ome5yMn1AyAJVyEzVRjdGxe9Xit6wAul3o588O57O7efHZPQAghhBDHkaCTFguoHcodcT8AFptJQSqFxTQxTIN4fk4m6OSHSzGsqqLjqlYrKdN5AFee2ioiN+bBb+6TcTpCCCHEMJOgk5YKqW0eHDE/msVAG70Ea24pZcn0nlcza8gJq6BTHqgmYc9T/547TR2gYz85hbkAOJMeLO46tjf2nNVzEEIIIcRAEnT6hNTO5Pa4H91mQvUizNmfZFRSVWVqK614HWpH8/yo2rFcT0aomTlLPb67Hk+5DwCr6UG3Bll3aM9ZPgkhhBBCHE2CDpA0klgjqtvJEevBYjOgah67Kz5AVUJVdLY0vEXp0lloZirzOEsqiC1/FNhywEyRW5h+OTUvmLChedNZPxchhBBC9JOgA/hjfnLiqhrjiKUrOpXzeGy/hiV9++GWzeRdfjGucHvmcRZrAjQNCscB4M1XFR9Tt2JPOqkPycwrIYQQYjhJ0EGtipwJOnE/FreDuKuYP21roiWuQkx3rAl3lR1HtCnzOKfXpv5RNAEAl6MTS1ItLFgQzCVqOSADkoUQQohhJEEH6Ip04U6ooGOP+dG9eazZ1053OMHB1Fx1H0uMeN0r9LhaMo8rGKUeQ6EKOnpvPY6UWo+nstuDbu9kXV3dWTwTIYQQQhxNgg7QHQhgM+yAqujohWU8teUIAPNnXgWA36JTt+kX7CvtDzpjpqtByX0VHToP4LCoCk5FbwUAL9e/cTZOQQghhBCDkKADdHWptXM0M4zFSJAqqOLF3W0AfHTueLyG6qLapIfZVtMfdMqmVat/FI5Xf3fsx+VU09QrY2rfrO2dW87GKZw7DAOe+jy8/N3hbokQQgghQQcg0K3G1VhSPQDUaYXEkwYTS3OZVuGl2FEGwOtOJ1057YSsvTg9NjwlOeoARRNBt0G4A5dbLTxYkO4KM2KvQGoEjdNp3w3bHoE134dEdLhbI4QQYoSToAMEe9QF2Z7sAaDZUFPNL5pQjKZpjCtTiwJucDkxdIM1817mk/cuwmJJv3x2N4y+EIAcazcArpj6XpstQs/vb4FU8iydzTDrPJD+hwnd9cPZEiGEEEKCDkAsoNbGcSZUF1aHqbqqyn0q8IwvHAtASFcv19KJ43C4bQMPMnE5ALlmPQDxuIXyZBJD09hzaBU8deuICDtmx36+W5DPL3zeo0KPGAkSra10/vphjKhU8oQQ5w4JOkAiaALgjPUA0GaoVZLL0kFnlGfUgPvPLBt7/EEmLgPAE98OQFTL5YKQWldnq9MFO56Ajb8+I+3d1bmLR3Y/gmmaZ+R4Z1Jt+3Ye9Xn4cUEeTa1bh7s54ixqu+8+2r7/fXr++MRwN0UIITIk6ABmevsHZ0xVdJoSFqC/olPtqR5w/8rcyuMPUjAWiiaR6+gEIG73MqNbVX1e9KYHK29aCWcgnHzn9e+w4s0VvNTw0js+1pl2xF+f+ffzbRuGryHirDJNk9AG9fOOy5IKQohziAQdQA+rqeWuyMCgU+odvKJT5aka/EATryLHosboxO1eKpOTATig92JYHNC2E5q3vuP21qfDxPrm9ad9jFAsyS9fqaWxJ/KO23O0I6H+WWnPhw+f0WOLc1eisYlUe4f6d1PT29xbCCHOHgk6gD3qBsAZUkGnR3egaVDiUUEnz5FHri29M7nFSaGzcPADTVqOW+8B1DYQJa65mIaNlBahbtIV6j5bfveO2hqIB+hN9ALwZsubp32cR99s4Dt/2c1/v7D3HbVngEg3R8z+8Rk7tQSHA+dn2Gnc281vvvEaDTs7h7spZ1VbuI3vvP4dmoJDCyuRLVsIuUvZN/5DhJs7TnrfZCKFvz38tsfsfvz3HP7s50gFAkNqixBCHG3EB514Ko4r7gHAEVa/UP2OXApzHNit6uXRNC1T1anMrUTTtMEPVrUAiysXmxEEQE94MaOq+rOuUFV32P4HSJx+FaU52Jz5d52/jrZw22kdZ9sRFer2tfae8mNMwzj5HTprOWJV3YBauovu+QPPnFb73k4yZWAY794Ypc3PH6K3M8q+N1vftec4F/121295fO/j/GTrT4b0uMiWLewf/2GOVL2HQ8lRJ73vK4/t43fffJ3mAz0nvI8RDrPjF39hY3MF3X/6y5DaIoQQRxvxQae9qwuLqS7O9ngAw+kiaHNlxuf06euuOmG3FYDFChOvxmWqEBFu9ePT1arJayO94BsFUT/sOf1f3I3BxgFfn25VZ2eTamNde+ikg5pjdXU03HorB65cxu4ZM9j6oWtOHHg6D3DEpl7LK2LqPs/V/e202ncyB9p6mf6t57n3TzvP+LEBoqEER/aoLkh/+5nt2jvXdR3Yxd1Ppmjc/vqQHhfYtoPuPPVej5BDKhgc8P1IPMWtv9nIlx/fyv5NKpw31/pPeLz2P73AzjEfo6liKfteqR/aSbwTyRihb19N6NvvHRGzJIUYCUZ80GltV10TSa0X3UwRKygGTcuMz+kzPk8NKB7jG3PyA171H7jTISncEWSsR63Bs8+/E2bfoO7zDgYlN4eaB3z9ZvPQg04olqSuI6T+HU/R3hs74X39Tz5JaO0rJA4fRksZOHbW0tN6aND7mp0HaExXdG6ylWI1TfYGG6jzn9nBqX/a1kw0YfDohsMEomd+Mcb67R2ZatFgXSwpw+TLf9jGN5/eMWhI7HniSfYtXkJk27Yz3ra3k4inOLCpjVce30fHkVOv1vUZ/cIOFu01WfpC8ylXC41QiJZWMHU1+D7myCPROLDra9XuVl7Y1corG5pIRNVyDoHOE4fIbc/XkrKq/0ctzUn8/uAJ73smJf/87zQ8Wk/DI7XEV/3irDynEOLdNeKDTmdnX/9/DwC9viKA4yo6N0y5ga/O/yq3TLvl5AfMKSKnugaAcG+C+WWz1NGTh/FPez/oVqh/Bbb/8bTa2+hvoDBgsmx/HDi9is6elsCAnFWbDj2DSR5Wa+G8sdBLwKVuO1y/Y9D7dnXuIaLraMDUgsksjKjxOmsOryHl9+N/5hmO3P456t+/jMT2NWq7iNPwem16ZlvS4IWdZ75r6eDm9sy/I70J4tGBn+zX7GvjD5uO8NvXD/HstuPHsvQ88QSp7m56/vgEWxq6OdR54tf3dAVjSRo6B4awTc/V8+t/foXn/2cHb718hL+u3MUXH93CLQ+/yZ4j7bDhV9B74tfLNE3cbSocTW0w2dK6+ZTaEtm+nc6CqZmvo448Ek0DK4+r96rQNC5dPQXYu2Xwilxwfz11+sTM1wHPaG77yq/45Su1b7+kQttuSMbfts3JlEFdR4iX9rRmqps0baXpTyvB0MDU6Pifh8BIve2xhBDnthEfdHq61CdFW1L9suvKyQf619Dp43P4+OTUT1LoOsFA5KPklnoBiGkuLsjJwYip8PRsx2ZSF/2zutNf/x8Emk90iBNq7trL5/9i8Ok/6syqM2gMNnKk98iQjrGzaeDgzrqTBJ1Urdqra2tBkG41Hpu2xsEHMB/pPghAqd2HvWgi86KqUrSvYQsHl11F01e/Ru/f1xDZd5jAd/8B/msCrB/aWJBoIsWWwz2Zr3/91u/44ktfJJx4+8Gtp8L/6ps07FBBR0cFnMAx45h+s76/ovXdv+4mFOsPQmYySXTXLgA6X3mVD/z0Na764VpePXDyAbpDEU2k+PBDr3HZf6/mhT0H+PorX+fJPU+x8S/1JOMG7nyHaveRIH/f0sTqve386edf5bev/BtNj3z4hBfvjkgHRT0qSOSHYN+WVafUnvDmLXQWTM98HXPkkWjuf28bhsmave1crb/JlZauzO3JwOBBd9Nv1pO0utHi6v9kKKec6Uf28Z2/7OYzv9mEP3yCKt7Op+Cni+CJT52wraZh8tqBDmbd+wKX/ddq/mnlRq578FU217bw1z/fygs93sx9u/ckSK5/B5MHzsF1roQYiUZ80AkFVNXBGVcX/xan2qOq7Jiuq6Fw56vSR9zuZWyih2RIfTr9/obvc133Ol6rnKbG6jx7x5B/GTb2HmZsi3rMVbvVBXZD09Cmme9sVOeqp8dUnzDoxILEu1X1JOCCcK76NN7VVH/8fU2TxnS3WmVuJRSMY3pMBZ3g9q2k/H50jwdHkXrLJaNWCHeQeOl7Q3oNtjT0EE8a5NgtgEmD8RQvH36Zlw6f2ppCoViSn7x8gMNdxwcjMx5n67/9DMPQcIVb8UYbAPC/9JvMfQ51hlizTwWhUq+D1kCMH7/UvwJ07OBBzPTKwJaWJkpCXUQTBp/63w28dqADI2Ww5tG97Fp3+lOw73tuD3taekkZJv+y5j/5U+2feGjVwyQTBk7dz8S8/6TbEkFD4yPlRVw9tZhQ4Wt8vzCfHyQbVWVnEIcDDZT29H/tf/XV/mrHSbRurSXmzEdD/Rzjdi/xo7qudjT5KQ4f4EHbTwkF8vsfaOYcV6GJhePsbVGTA5qMOmx6CDSdhSSxW3Re3N3KtQ+uoy0wyOrLrz+k/t79LOx/8bjB6s/8cAu/+9f1/GLti1D9H7hL/0aRx0rSSPD9v32Sr9nDTKtX9w06QU9p1P38/hNWHuNJg0Sq/3tGyuCNZ2t59Ym9vPzre/jTA1NJ7nv+uMcl4ik2/KWOA5vazslFP4XINkMKOitWrGD+/Pl4PB5KSkq4/vrr2bt34Kd70zT51re+RUVFBS6Xi0svvZSdO99+0OgTTzzB1KlTcTgcTJ06laeeempoZ3Kagr1qnIA7qi72h2zpoOM7/aCT41Xr8sTsXnI6WskJXUus4zLcVg+Hehv4Z3eKpMUBB1bB3qEN1vX3dJGb/h0/tU4H0+SN3X847n6ppMFLv93NE9/fSCI+8BP8jvTFa9FYVZ2qbT9B0DmwCn9cvUVybSmK8wsACLU2Hn/fYCtH0hWQqrxxUDieqXHVhWA/oqoZ7gXz8Vap6sh+96UA2BIBjNCpT+Hu67a6fEopU6pAs6rAsrFl4yk9/j+f38t/Pr+Xrz+1/bjvJZqbafOqLpiS9q04/Oq5/Ds3w56/AvB/bzRgmnDJxGL+4/oZAPxqXS217aoyGN0+8LhXGL9hzsS3iCYM/ul/N/DymgZ2rGlkzWN7iQTfvovlWOv2d/Dwq/UA5Ob0EHaoQcMVATUQuNF7gLnR1Sx3qa7RuRYH918Y4AWPzsQjJlstDlIv3guB44NW04GNOI7qpSupD/DBn71E62ChIs00DBqbVWKuHGVDwwBNJ9jYX7l5eU87s/SDHI5PBmzY4ukKmeYkFh7YLVj33BYSFjeOaCeds3/Jvjz1eoYNH09eN4qqfBcNXWG+9uT2gSGhdRccfiPz5ZFHv8iV/7WKrpB6jYPdMY7s6SbQESV/Xw+6rQdLwRqmzHmE0rG/YGduK2NaTPJCELHY+fXFaiPfwFsxkpsfO+68/eEES+97iff+6BU6gyrQv/7sQTb+tZ6tqxrZ9eaV7N3/PX765IPQeTDzuHAgztP3b+HNP9Xx/P/s4Pn/2Uk0OII2/RViGAwp6KxZs4bbb7+d119/nVWrVpFMJlm2bBmhUP+F8vvf/z73338/Dz74IBs2bKCsrIwrr7yS3t4TD4xcv349H/vYx7jxxhvZtm0bN954Ix/96Ed54403TviYMyXcqy6U7rD6uxbVP/NOgo7bq7oO4nYvicOHmVxaRLz9Klp3/DNmykVvMsz92jwA2neuPvW2JsI4uxPUV1/F6/O/SUwbT3kX1HbvG3C/VNLguV/sYPerzbTUBmir7++qiieNzJTy982sAKCu4wQDPXc9Q8wsZsfUTzEnmUeOru4X7xhkkOpRM64qvaOgYAxeA2oSCao61QXJUVGI1aEuPIcO+2kyVXBqqh18zM9g1qeDzqKxhVwwvn8w66bWTW/72K5QnMc2qCrNqwc6aPEPvIDHDjXQWagGjycmvI4roio3PYlyzD/cQvLF7/DsBlW9uXFRDZdPKeGSicUkUiY/Xa0uZpHt6lySFnXONS2N1OuPcNlEN9GEwRPPqfsZSZO9r7cwFIFgkEOPfonbLM/wxTkWZs14E00zMEKTuSSgQtfsrQd4YXMNFSH1ehzZ08WT63/BFW9qfOe3Kd67VqPJjMLfvnrc8bv3DAyLUxpMLNa9/N/rgw8+B4gfOkR7jhqoX2P9Ow5T/Xx62/rfU6v3tTFFO8RrybkAFHbtzISdA4faBxyvu069t/RkA+vyLRzKUwPZ/b5xVOzdyq9vmY/dqvPSnjYe33DUGk2bVgLQWLiYdtNHldHIFf4n+eMmdZ/W+v7K1PieidR0TSfHlsPW9s2EHQ3kGAaffUu9f3fVLGKT6zO05oE7qvGf//sdHtjwn+zs3JkJV89sa6StN8b+tiCf/s1G9m1tY8vz6rk6HTuI60FcSQ8cuo2N//rPdP7Pz9j7vZ/z+3tW0VYfwGrG0TA4uLmN333lRXY/thoj/cEgEUvRUuundms7u15tovlAT6Y6ZZomwe4Yh3Z2svmFQ2z4Sx2H93SRiPV/mDFNk2gwQU9rmFjk+JljqZRBb1eUYHcU8wRLNBip0xs/J8S5yPr2d+n33HPPDfj64YcfpqSkhE2bNnHxxRdjmiY//OEP+cY3vsEHP/hBAP73f/+X0tJSHnnkET772c8Oetwf/vCHXHnlldxzzz0A3HPPPaxZs4Yf/vCHPProo6dzXqcsHkqQAzjTQafeqsrm76jryqcqOnG7l8SRHXz4Y9ezuaGbaMJBMjwam2c3G5wuCEJ37SaKT/G4Tf5DlPaYtJQtJOwuZeusO7ho/+Osmbce4mGwu0mlVMipf6t/TEjY31852N/WSyJl4nFauXiiGjvU0BUmmTKwWo7KvfEw5p4XaC++hbaSC3CFkvhSDxLHg9blJ2WksOiW/vt3HszMuNp72M7vEx181DeKabEQFZ3qtbQX2LE61S9QXzRIvVFGhaWLxtodVM289G3PP5pIsbWhB4BFYwtoq++C9Nha14469v3yvRTfeDP5H//YoI9f+Vo90YR6fsOEp7Y08vlLx2W+H6htJGUpAdPg+7M6+X8H1Wt4KDwGLRXDuu4/+aNZxO9yP8pl4y9H0zTuvGICa/a18+zWJr569WQib6mZVn+fqXPVFpPph0ySmNxY/SJ1nVeQV58E1Ou2a10Tsy4fNWBdplQwSLz+EK7p045r/47nfsUNxrNgg4MHnmBlZQVoYGtdSrynFDQob9tPbjhOz17wLGmk117Jxv1O/vFVdd4L9plsmO1k1O5noacB8vq3NwkdqiNhdbNv4iVY427QbFwV7OL/3mjgtsvG47RZjmtTcM9BAumZiI/b/0qucxR5sWKC6fdcVyjO1sM9fM3WwJvBq9EBT+gAIXc5CbuHjXsOMH1qReZ4ajp/DkGbCkAt3lr1s/HUEFj3DBNv/CT/vGwi3/3rHr79511cOL6IUbnANlV1+VrzJRQzi/vtP+OL1if55BvL+cxFY2mtU2E/aUliTVm5vPYjXDP6Ff5V30q0t4n72nvY1zCThLWb0OiP8MlenefHL+CmjW8ybh98e9dv+J9dv2F0bhXv84xn1e6rM21u2tfKqm29oNmpbFzDe/b/HkOz8tb0W+kqnMbG2KfZ+/fD9KSn3zsjHcx66yekLA52TbmFcE4ZL6022Pn49zDznXQ5Z5E0Bm4a7HBb8OVb6WqLkkwMto6XocaUmSYmVkyt/2dl1WNYtCiplIZpWEjpOZnvaWYSmx4EM4FpJDGxk9JyMTUHuhnFqofQiIGhjm2YblJaDhqgaxF0PQZmAowEmFYMcjCxoWlxdEscSKAZcTANTNOFgQtNA01PoekpMKJoZgzTsGLiAqygG2gWA4iDEVNdh6YDEwdgousGaEkgjkkCDBuYDsCCpvU9NoFpJtTrYdrBtKNpJppugJbCJA5mEkwbGHZAR7OYaJqB+h+rHothB9Ny1GMNDM3AxATDgmba0ABNN0EzgCQGSTA1NMMKpgV0E10zMDUTQzMxTRPN1MG0HvVYUz2vmch8D1NH09XzmpgYGJhoaIYl3SZUm1FtNswUWroZGhpoGug6JqD+95topolmmGi6hqaDqanvmZjpNqt149A0NE1Tz2saqls6/bbTUIc2MTFNMEk/X/93j19rTjPT3zPT9zLVv9MfHKLxd3cj4CEFnWP5/epTUkGB+mReV1dHS0sLy5Yty9zH4XBwySWX8Nprr50w6Kxfv54vfelLA2676qqr+OEPf3jC547FYsRi/dOiA6e5eqoRU78Q7PEwaBqdLh8ep5Ucx+m/NO5011XC7iHacIQPXlDFB+ZUEk8Z/Hp7Cz99azf2GivshKLQXuKJFPZBLiLHamp8g9JuiNtVGDN1KwXmDZS1OTAbN6ONWcqO11uof6uDlAZBq4kvodHT1V/56BuIPK0sl4qW1RRao3QmnTT2RKgp7P8FePD1Z6gJRYjZ8wCIJWbis5m0A56gQVOoaeDWGJ0HMosFPrspxrPRt7hu4himdW+mslNVuByeOJpTffIsjPeiFU2E7l0EGwdWpAAO7+midks7F354PNb0a7O5oZt4yqDU62BMUQ7Nu9VFcO5+gy89ZZBK1dH9+8f7g45pQuMmKJlCyHTwv6/VA3DZpGJe3tvOk5uP8LlLxmb+UwYb2oESdLOXoFVj56hOioBoqpDPx+/kG7b/o0rr4GvJn8KPn4IrvsUFsz7O7FF5bD3cw6Pr9nPF3n1owF/n61y5Q6MgmKS8CzZFnuSB625j9f17VNt06G4J01IboHycL91ck8Of+xyRjZuo/vWvyFmyZMBrYtStU+20l/CzPANTg8sSFhYlm+jWZmCL97JuUif2lMZlOzUKmt+it6aSKXWzsCfVPlQFQVgVyeODBKFtz4Cgk2zp5HDV5bSWvDdz27iAwXPEeHZbEx+dd/xCgB27jmBq1dhSPfy5IMrylh7yYhCJ6ZjxOK/sb8c0TfIs3ehGDWhQOTefnv2d9FJDU23DgOP1pANSwNnJZ2JO1tk7ietB7OTStquVmlSKTy0dy4u72nizvosvPLKZPyyupysZ5O/uMt7IbeWGae8nVf8aOa1vMb/nb7xZtzQTdN6s+DPT2i7EFyume12Exzyb1YVg7m0Yj/yZIyWLSRo6FqCj8HrgTaY3GFzTGePvhbnUB4/wYPAIH2AVduut3PKxj9Ow4s/EvOPx9DYw6tCTBPNySIVhwp6H2T7zM4Q9k+jJm4hmpqhiN2ONX/LnpTFaLFZc8e9TGb4WzXkZraVL0z9osMf8OKOdWIwovbk1xMI5tIVTgIZmpnCF28gNNaGZBj2+8cSc+RjYMxcjAEsySsrqJGk4SOLIvO8ANCMJaJi6lbiZN+B7mfeb5iRupj/waQw4NkAKd991a/C+gb7vDfLYvvPMfP/YX39HF5Tert/h2O8fXaQa7LmP/b7lBN870bGPffxgtx/73IN9v8+Jimd9ueBEE/9M4GRLPZ3sscZJnncYROJnfmbq0U77am6aJnfffTdLly5l+nQ146KlRZXiS0tLB9y3tLSUQ4dOXP5uaWkZ9DF9xxvMihUruPfee0+3+QCkjBRaXIUSWzJEqqCIlG55R9UcAGeOTSVeE4JNXSrBaxoOq4WLRy3ip29BbaKBGDoF9PLqjl1cOGcGf9rWxCNvNHDh+EKumVnBmKKcAcdtat1KaY+FZL7asqKiaR1NFUuZ0rqUYMM6PGOWsmGLes222ZIYGszDysGGAPPTx9iVDjofd7zE43/9Gf8vdyFf7/kUtR2hTNCJJlLs+vvvqIrpxB3qImyknPQ4xwHt+MIm9f76AUEn0VVHi1X9xjAT+Zgm7ElVMb33DVzpXgy7o4v/KypiXN5E8vwHsFROgu4/YempHXCepmny0m92E+yKUTbWx6SFarzE67Vq3MeisYVomsa+7n0s3GNw19MG6Z4i4gdrMVMpNIsF1v4XvPwdmPJ+Hq24F38kwZiiHH74sTks+O6L7G8Lsr3Rz8yqPAB6W3oAiFjUa1S47D3wAhi2PFrtlfx2zh+4wfIi1Xsfht4meOpzUDmXf1o6hi8+uoVXn1/PlYaB3w0d+TruObOJvrmB6YdMXp9m8PFtf8bCBAKawSGrwYy4lV3rGjNBJ7hmDZGNqsvJ/5e/DAg6kXiKmuBW0KD18u+xet93wYjz2bYjtO7vpbsACuINxL50Kz/b/TN6JpVwxd+2cahmOYZjKgmrHa2iAGtDC3qLDqXQXvcWxRP7P5Q4u2ME0xvWFtmaCQRdxB15VKVSPPxqPR+ZW3XcJ7Wehi6gmoiuql897h5q/BCz59Fy8BA/W9NCJR3sTZWDpuOMdDDusx9m771/Uu+btoGDnYMxO9jB6mjh9mUPMHrvY7yyr47RPTPodlQQ3bMH17Rp/PdHZ3Htg+vYdqSbOzf8N+uqKynvNPn2Cy9Qrx/Bcumt8OwX+AfLS/z49XomHVJdZZ/983bWz/YDN3PIciULva+g5RYTti2F5NO0VV+Uacscbwl13jLGBFq4560Q3xzTyR89ufxXYT7PeS381f8dOlceZL/3Q2Cm2FT2G55fZLK+8d+wWa184bLx/GrbQyzxt5OwhJnpeg4brXy6qICQbseR0rim28LSwNPUJXfQ1HM9iXgnNYdfYVL93szgbkPTCXjH4M/NxzBbiOotBF0pWt0QdoBumHhDHgzsxHUbSS1FQvOT1BN4gw68iQIsuIhaDWJ2k6QtRNwWxJoyyIl4cCTyMLCRstgwLCkStl4MSwhHzIkj5kHDSUqzYugWkrYQMVsvlpSJM+7ClnJjYMXUHKQsCRLWEIYexZ6wY0+4wLRjalZM3UJKj5CwhtENA1vChm7YMDUnpubA0JOkLGHQ4liTFqxJGyZ20B2YaBiWGIYeRTNBN6xg2tCwgWkHPUFKi2DqKTRDRzfT3zdV8DO1GIYeQzN1NNOC1vd9rBhaAkOPo5EC05J+rBUNK6ZmYmpxDC2JZmpgWtFNHQ0LoGNqSUwtgapMWNK3q+ObmpH+fjJdodHR0NP30zC1BKaWTil9xzStaFgwtZT6PqZ6XtLtRofMY1WCUcdW7dGwqJUR9HTVxDRAFYvSeUtLr5zQn8A0U9ViNFNLn2/fY/sqLlqmWmMOSIx95Z2+Y5lg9iUnc0DtJnN/Lf13JgAf9dhztaLzhS98gbfeeot169Yd971jfxn2XehPZqiPueeee7j77rszXwcCAUaNOvnS88fqiHTgSKrQYEuEiZWrTqR3Mj4HQNM13B4boUCCSChFZOtW3HPmADCpYBI5thyCiSAbc6u5MFjPrs2vMnbseL76xFuE4ynW13byXy/s45Ylo/nW+/u7MJq6D1DRm0tnPpimwdjIJppYii9aQnvD03iA3pYwdqC8OodYxIDGBN2dR1d01IVlf+olflNUwLLQXuhRKyRfNkndZ/vhLi5hE+G4hZjdl3lss3UWRbyILwS1gXouov+i0BI4hOHQsJgWzKSqOD3dWcUXOw1agO4ceK0lSEvbd+mdnc+0nb+i0qsqgYUxtfCf16nK9a11AYJdqlrX3dKf9F/Zr7ozFo0tJO4/Qn33Ae5YrULO2mkai/eCLRYjceQIduMQrP6ueuDuZ3l5/8VABbdePBaf28ayaWX8aVsTT25uzASdUHsACiFs8+OKLuT2m+7ioeeex9TtLO9ZzcUXXUpF3t28MXsJL750D6M7arnhlftZ/v6fUO5zMi491uhgucbsohl4Fi/OBJ1VF9g5uGkXMAF7uZu3An5mxK3s3dDK0o9OxO600P7AA5lzDb68uj+wAZu3b+dCrY0UOu1lBUT3xBkb91KwPsw2l3qPjFk2hytnzeC3+/+Xx8a3oc80KI20E3EV03jVHUwfHyH4o58y9nCKnnKduj2bKb5KPV9tRyf5fpNAmfrAMfuiKnY88hotZQuZpnXzdLOFN+q6MgPY+3S39kIedDs7yE16SNnUzyvqyONffvYiu13VXO9qpKOnBABbyo9j1oexaSsBcIT6/48bhkncqn4WrtwwltGLuTq/mj9v/B70zKDHN57Ixo24pk1jVIGbh26Yy9efuYt1riQWw+C2P1uY1JTC/ZfdJG77Ffrf7mFMopXI9tdIxqdhSUYo6m7lqlcjvLoE2gP5xL61FYfbRu+//AthVzE9zsp0WR7M9hgbquYzZtefCHSUMGqunX+4+Gvc/8YviTgaeSrlwdU0C3wQsaznpZntlDVdgtVi5Rc3zeOSicVcNuVr3PSX20jY96KGjavu4lR4NBcVf5m8UWU8EYrztx3NdKTieHLcTF8+ho5oD1pvL5ZQjEqHgzk2P1MtFvJmfATnmCWYuio1uK1uXHoJ+xp7GFOaR1X6Q9DRr+mhzhCGaTK2OPe436kpw6QjGKMo14FFP/73bSyZQtc0bJbByyqGYaIP8jghhiIQCPCF33/7XTv+aU0vv+OOO3j22Wd5+eWXqarq3xKhrEx98j62EtPW1nZcxeZoZWVlQ36Mw+HA6/UO+DNUzaHmTNCxJkP0etUv8Xda0QFw56ljxO1e/E8+mbndqluZXTIbgB3F5QCEG7bwrWd3Eo6nmFzmYXH6YvKX7QPX2WkKtZAXViEiSoq9ZTXY4wE0dBoPd4FpYutRF5qbYz9gbLFK2OGA6g4wDDNd0TF5FTXos96axEdwwBTzA/t24tUitMXtxBz9QafRsgAAXwjquweudnw4rAbLWBJe+lL8X7qr0f3qQr1/zKVs3XUDOQk1vTjgHc3PDv6VBDBaa2Xroe7+59/YP9i5s1m1a29LL1saerDoGu+ptlK3chla0qQ0/bDfXq7TWKA+KsS2b4I/fkp9wrCpKtUNsd8ztiiHD16gKhYfSv/9zNZG4kkD0zQJhdXjA84A1476DJqmkZOrbhu/KU77Vdfz6pKZ/PfDn+Mxs4fvFRawb/cT2AIN3LR4NJf0qPE5B8thcc2luBcuAmD2IQ1r0uRgVI3RuOrSai5ZXEmHbmAmTf771xv47U/vJLZrN3GHhYTbRqqri8jWrZnX4ci2vwPQ7JrI6x3bWLzL4NsPBug64MTvGwvA6Mtn43P4+PDEDwPw6MUaVufL6uccGYdl2sWAWgxwq82BpXM/bb3qk9QfN7xCSY9OxKUuxGWLZ5Dbq6qwY9O/Jh74+/4BM53MVIpQVH1Wavd0ct8VP2BihRpvE3PkYW9rZExRDv+2ACIx9Z622eJouo7bpgK3M+kmmlCfToNN3Zi6Fc1IUVip3nd2XyWzClTA7SqYwta/P51ZM8l07qWzQL1GX3zFy6QmNXtpdLPBzrY96LM/DsAlCdU95u09RE+OicWWwJV+vzYd8GMmEgRXvUhriap7Vk0poHJiehr8WFVV8x+2Ynz2NV7QlhBqew+6YeLbvBC/bzy6EeePs5/Hi4U73/dVnvj8Ei6ZqD40Ta8o4aUbHmGscRtGXAX7VM9FPHDJz/nBhy/mS1dO5NvXT+eNr1/BX764lHX/73oe/eD9rPrEy7zw2Y387e7t/PL2jdx+6ytc9qlnmLPgVqYUT2dq4VSmFk5ltG80pR43F02uOC7kAOi6xpjiXMaVeAb94GjR1Srwg4UcAIfVcsKQ03d8Ic51Qwo6pmnyhS98gSeffJKXXnqJMWMGbocwZswYysrKWLWqf6GxeDzOmjVrWHLMeIOjLV68eMBjAF544YWTPuZMaOxpxmqmu64SYTpy1S+3Y1dFPh1943Tidi+Bv/4NI9y/bsu8UjXjaleOukiMM+p5fmcrugb//dFZPPgJVf1p740RS6Y7WVMJ2qK9uNMbkAZ1eEqvJDeoFgtsC5cQObwfLamqIjWxV7ii5QfqoSFVIj3UFSYUTzHatYeDNvWjr7fZmKXvHxB0OuveAuBwqgzDku7b16AzUUnM7sNiQktr/9oxJCL0HI7w779N8p43NXwWk2tmltNGPv5wHt15EwgWfQTdtOB3qC6OUE45rc31/Dg/D48WYfcBNRvJNEz2buxfvbe2tgeAR95QF90rppRQ2vwyexM9lHVD1FVC/ZiruLjx49SNu52OwunE/vIjCLVB8RTevOR/AXiv5U1+fJkNR7p77aIJxRTlOugOJ9hQ30Wqs5O4rl7bkDPMdTNVeat4rArvSXsxhb1Q0mPygQ061R41tuW33hxY90M+MSOHmh7V7gPlGgvKF+CaNRNrSQmuiMG8/RYiMTU7KfXdL/LpX3+dUUa9OskdOyl/RL3//zTXYP2YJL05lWx8YicdR3oxTRP7EVUPMKqX8Hrja3z6eQNbLEls9nswLA6cuTYKylWou3327Xxmxmd48IqfctPPf8bY2cUYhslr6xNEcp244lDf42Sc1shvX6unOxRn2841uJOFmLoVq03DW5mPFfXesobzses6rx3sZN1RCx9ufvp+ok4VjEpHVXBxzUJmT1JhOObI4xpbPU9+fgn5vXtJJtX/LadbXRhLPCqUWM18djSpYzbuUF2Y9lgns8bOyzzPJ5dcT8LShKHbCLfkc8ljF/PeJ9/L3S/fhanBJ1qCzN3cP17PasC+1/6KNu+f1M8wqf4vegN1bJpgo+S3j5Lfs1+9xzcfIfTGmyT9floqVDCdtLCMiQvUh6xiSyEdTi+WWITAa+v50Yv7SQancNM6L2HftQDsKllLyBHg2gkf4H0zqplR1f/hACDPbefJmz/HP9b8lOrQCn57/Xe44qgB2KACx7QKHz73wEHIQoh3bkhB5/bbb+d3v/sdjzzyCB6Ph5aWFlpaWohEVNeIpmncddddfPe73+Wpp55ix44d3HLLLbjdbj7xiU9kjnPTTTdlZlgB3Hnnnbzwwgvcd9997Nmzh/vuu48XX3yRu+6668yc5Qm0dKQrJmYKSypKiyMPgNIzGHSSJdUYoRCBF17IfK8v6GyOdWAAUzR1Ef/kohqmVfgoyLHjTAeR5p5032X7HqIhnYRdVa56LTq13vJM0PHHa2hbrdYe8lhasBRWUJBQ4cEaT2GaZqbbqqZ4Q6YtcV1jnH1XJuiYponZrgbMxm2q6mEaEVI+9Qu4s0Ttwt5z9KKB/iM49jqZfAT+6ZU2fv7CCm5Mqu+3dudQO0ZdEHaXrKer8mFABZ28EPzB68EE2g6ptZYaD/QQDcQx0/29RiDB/pYAT25uzLxG9DSwz26jstNk+/Rbqa15P5WdSzAdU9kz8R+I1B4Gi52e9/0Pt6/ReJy5bHytGPNznyOeHitm0TUum6Q+db+0p4344cOZ6lXcZmFGuqLgK1fhx/rBGyn50f0AzKnX+O7cbwLwl9wcOrY9Qu4DszEDqtrRNMrF1MKpaBYLvuuuA+DCfaPQTDsOI4CtfjfxgweZ+aZadTcvVU1Zj5OoXWd79QcIlv4rG+Z/ne1dVTz/ix3sbg4wPaleH+eURaS278ITBc3rIfqhOwCompSPlv50nWPL4YsXfJGLq1QF56KPTcDmtNBaF6DpAvWzSLQ4yNNC/Pn17bz/J+vwdu8j7FYXd1+JG03XsOaH0IwkpuHmxlnqvfD95/ZiGCaxaA9PvLkyUwH61BW3ADB33GwAYg4fY5P15OfYCbdsx5LMU+daoJZvKClTi2pquodXD6qVpHdsU2ty6akOpo1dmHl7+WZ8lKWj1Lo//vz5FLRHOdx7mFAqwpRYnE9urcIejeCYMoXuxer9GdjwBpROxRy1iNaEWrDTF6iHmVMpnTqR0gL1AeLItiZ6n3+OgHcMEUchVoeFsbOLGXdBCRarjulP8Gb1hQCsW/kk+9uCXOqvZc6BCYRyK0lqYV4b8yKXVF3C5+fexYlYdI1/XjaNv9z2PubWFJzwfkKIM29IQeehhx7C7/dz6aWXUl5envnz+OOPZ+7zla98hbvuuovbbruNefPm0djYyAsvvIDH48ncp6GhgeajlohfsmQJjz32GA8//DAzZ85k5cqVPP744yxcuJB3U2u7unhaU2E04FB6avmZrOgwUa1v4n+iv/tqWuE0nBYnPckQB202xmgtVLlT/L8rVRVB0zRqci1opkFTjwqRsbcex9GrZ2ZcxSzQ5s7HEVEXgFhkNN3bVSUmZe1F/8zfcRWoT/gOU6e1O5qZcdXqHDj4N99RS2NPhEg8xZHuCBUJFQYShvqFbJhB3kqlx8yk15lJdXUSjKdHGfc0YOnV0+3S8fV24fv+vzHFEqHTnILfN44UcTZX/YULI+oxMUc+RSE3QV2j0Woh3rofwzB55ln1SbusbSOakcSKxmd/8Sa9sSSjC91cOK4Ieg6xz26jqquQUE45mpmiMP8PJLQAcUceDfp8AjNu5ronagnmP8B/l7WR02DD0pvk4Of+ESO9cvF7JqtxIy/vaSOx763MeCR3fmWmzO8rVhfkMLkULLsa++jRmPE4o3d0Mqt4FglN47FcJ7F2NSOmwwMTxy7Alt7g0veB6wEoDqlqjrf7IBpQcIEd97R8HNEWTN1K7dQpjFkS4MrW+TiNEjQjgWam6GmL8Nqru5mgq/fqFpeVOQdUCMy98CL2b1LdOuPnlpzgnQi5+U7mX6Oqrz1eVbUoOaIT06AkeojDXRHKY52EXSro5JWm3zfjRpMbUs/7vnIfOXYL2xv9PPxaPf/++D9h6bERc6hKTUGpel/mF3jUtFPdRntHBFp3sT/UiM1Qr21+qbp/zuhKLEn13t7ToN6PnY1qDR6ddvSSSf0noOtM/MdbVfvzJvDfu3P5TbKAH7d28vMjbYR29Kj3zDe/SeFiFe48e46QSCWIz72T7pQau+cN1FO+SIWW0VdfAEB32E7XS+toGHU5AONmF2NzWHC4rIyeobrbEqNVZbnqzZf4t9d/zZc2/Z660dcAEL+gnd9/6FEevPxBfI6BlRwhxLlhyF1Xg/255ZZbMvfRNI1vfetbNDc3E41GWbNmTWZWVp/Vq1ezcuXKAbd9+MMfZs+ePcTjcXbv3p1Zh+fd1NOtZvHYEqqacVBTnzaP3bn8dPStpZMqqQZNI7xhA9EDB9ixtpHmvQFmFavNPt/yFaFrJn/4QH/ZOt7QwPce+Rpf3fh/NPZEIBGlefujA6aW5/jsvGdKGUFdnYORGEVXUm0m2uwpQnMX4BozBx3VRbCnvpudTQE0ezuNlihW02SpQ3XL6M52wOT5nS1sbuhmgqYubvGkej2SWoSdphrn0+mdgomGLwSHAumZdP7DuHtVMPjG5deSmjYTkklub3ud1rIrAXDm/J1/8dcRTRTi0NWxKpPq4r/Hbqc81cTS771E6KAKY2WtGzIL9hl+dQ6fWFitxgR0H2Kv3U5JSH16L3SHuN75KA15amB8Xd7lfGjbQlosT2PNqWNs21HzOuuaaf3OdwBYOqEIm0WjtiNE16bXMhWdCyZOydzdV6SCjr89gqZpeNJLJ/S+sIqbp90MwOMFxQTGq6UT6ks1Fpb3B3TH2LG4Zs2iu0AFxLyu/TjyDXInHuLLl+zmtXFqXI+9+gIOVt5MzPTitDbhbfwaBV2qstb8pqp4dOeOZ33nduYcVOcTmXUZwe4YNqeFmukn34Nt0sIy0MAfdBKz+xjXBAdNO+P1Ri6aUERhNErYrcJSfpka65E3YwHe9M84XNfIZy5WY4Ee+vM61if3UNWtQotVN3DmqPeuxapjtauu0lgiH/7vI+yzW9E09dp6qlR4to2ZhDOq3rudHemfc7r31G7pAM/Arh1vaR7FzgBoOh3+2cw5vJVLwyGs2gLMRBJrZSXm2Cm4Jl6GoVkZfyTFzrbttNlVsHNGOohag0yecSkAZe+/Ene0HTSdg3lLaC+eg6bDnGX90+0nLVL/Pzx6EYc8ZTiMJItadtHqmUnUVUSOz85dN9/EuLxxCCHOXSN6ryt/UH2yt8XV+Jl9Zg42izZgPZnT1bc6cjSuk7NUrY+x5q6fseaRvfztZ9upcapP2M1e9Sm6PNo/5qXn97/HGYuwuHknTR0B2PU0e1MhSntM4jbVdeX02PnMxWNpdsbRjQTgpD6musSiZeoTrFY+HbeuRuvWHQmwq8lPrkdNX14YjnJZeAbOmEmTzaRaa+MXa2vZUt/JeE1ViVKJ9EBtp0GLxcTQIGlxE3EW4gtBXUANSI62HMATVkGn1TOFqjtuwwQcdUGCnlFYklE+6X6U94bCHNEryfOkFw1Mqe6QPXY7o7UW6IqSY+pYEyHyu/fiTg8YzTc07FadD89V59XhP0SXxYIDtV1D1dhcmi1jKCxchWak8PvGY+8IYPdtBeBr3o8A0JQPBibdf3ySwF//isdpY8EYdeFt3rOTeLpbcPmc2ZmfhTdd0Ql0qFVk+4JOcO1aLi1aTGVuJT1mnKdef5rm0gW0lL+P4l1T2L76CMn0IFv9qg/h940D06CkfSvtn76R/5o4lwabDb9PbZraEJ/L1h41BWppzu/prI5S1KEqdIUR1Qb3hIvYvWcdo9vUFNEjqIvyuNnFWO0nX4fJ7bVTOlqd35FRs7GY0Bhy8M2FFn7zTwuwdMczXVd5pernXjZ9Ad70gOSm2i4+fdFYxhbl8GHfI7RbLRSFVLeVxztwgTBvgWqvPZHPod4W9tnspCwq6HhHq+ewTZiZCTo2v9pFXDdUiM91BUA//lfThDkqzB1yXwQffwT+4XF6U/PYP+4DvDTxK6z82mv86bEudk69AXcMdm14jtptKkT5ArUcrLIwqUCFY93tpiRPBejD6WrO9EuqKKzMzTxf9fRCXB4bRA3uv/jLfOmyu0jdeicNE1V35Pz3jcH2Nq+7EGL4jdig07lyJXPXq190tmSIpMNF0OZiybgict/BYoF9+rquQv44Zf/6TcKzr+BgpbqQJeMGhY3q03GbKz1brEXt6WMmk/ifeRYAu5Eksm8/bPglf/B6KO3pr+h48xzMrcmnOb+InKAKJhEjD4Cy6vQxS6eTY1FBZ3dtDx3BKE7fZgA+usvFjO//jX9YY1Brt7LAVsuu5gDrt2zDrcVI6TZSKXXByy90YGrQoadncblLyQuZPLVdXaQP1e0G1EaIk8dOwpg0l7eWfJV94z8KQGXHy7gtah2TVME48gvU66un1BiZvXYbC33drBijZt55ehs4XJTCHVGzryY5Hdx5+QQKcuyQjLEv1ome0om51diL0XMrOVJxFdfHm3D3qvO7SE+AlmRm8UyKDvkxgfWzprPq4jt4+dIH2fxH1fbLJqkqhtGbAE3HxKCooL+b1VPgwGrXSSUNulvDOKdNxVZZiRmJEHltPXddcBcuq4ui0DR2T7kZt+Vq6v4eZO1j+3jjGdUl0+hRFc2Crt105fi52/EX/phQ5/bl8VPwFNhJJkyioQS+YicTfNu4Oq+Tos6+90QpoVQ+LZVTKX9LdfnaZ86ibqcaczVh/olnJx5tzCwVTDpL1GD37qAdZ88B4kfeJMdvZoJOX0XHO2Ey3kA9AB3tNtxWnb/fVEbIsxfNMMlJqJ9fXmnugOfJK80DIGnP46VULs3xXAyrqpL6xqlwq5dOwJFUXVVlsRRfemI1pkUFmerSwcPD5GtmoRlJgs5yth8YjTnxKt7aY+HwqCtImtbMkh3tRfMJuUvpen0ze15Tr1dF82sEJpRhs/QP9q25cELm3w4HLHjfwMkVFovOxPQaTjdVlPDtr36EQM0VxEw7eaVupiwpP/kLLoQ4J4zYoNP+4INUt6lPnrZEmM6cfNA0lk07tYvG28nJU0En0B7h2Udaeav6Y5iaBXusBwDrblVJaEvPAqJlOxv/Ws+v/t9aOmP9FaXSPa9woHUrb7iclHWTGYxcUOhS0z5HjyE3dCRz/wQm48akxwqU9ld0Dh/uxJa/nqgjQK5hMKqpkMNVlzGpKY9am40PlqYvCOnxOUbBeAxTXcCKy3wsGltA11FBxxeC9Y1b2HSomyNH1PN3eK0srSrgyf/cTKe9Gs1IUHPob7hD/XuW5VROoaBCnV8cdWHb47BTGGnAsk+tQ+OOtPHQNRYs6TCwMNfG7Zepbi78RzhiszK+bQyGxYktHqD8grGMuupODqaWkky9BkBudAzORC63TLuFyI6dbJ9+K/l8HpuuuqUaY8UQD/OeySW4iWImVLgxrFH0o6bT6hadkhr1mrfW+VX31ZWqO673+Re4eszVrL7mJXpK1GrMJTXWTPB4a/UR/O0R9m1WlYvyltdZe2UZgZQap3TT1JtY9N4fMfaC/vfcvGvGoM/+KJMcMeLO3kzQqLMs43tdG7mgr9tq7nIivQlcHhtVk/M5FWNmqmAScowlaXGQCFihYz/73/gxhb3OTEUrr0QFHT0nh6S9E0syimHY6G4OEVv1r6zKcVHih5hd/fzyqgcOrs3NS1czHXnsiecRCaqAZUlFsXvU/zmcPryooJMX9aH5fodhVe+3cRMHXw/LVeRlTEoNWF7751ae++EbHCy6FIClHxzN5x64VIU5TaO++mpsDaNJJgxywofJ69mPc9asAceruaL/64XXj890vx1tymIVZnpre/G2xNiySk1VX/yBcQPeJ0KIc9eI/J9qxOOYwRBJa/8aOg3pLqErp5yZoOMtcjHtogp0i0ZrXYBwIE5+mZuFjb9FM5Ik2m0UhippNVVVKdW8k60vHiIWgz2TbiBlV5+AZzSs41FvLpgmZX4tU9EpK1FhoWz6DDzB/t3EOy0m48s8RHft4o1P30xPixoPU2Z24Ch+HoAvdfRQH5nJ/vEfJuq9BiOpUWE9gK7BRE2FFq10CoauXpOi8iI+Pr+arvSKm31BR3ce5o5HN9HZrsbV9PhymG9RO1LnlbpY2vgw4+r+TK1HfYpPmBYqR0+mcLS6MIet6kLZYrXSYybojqjqSsjaRnRsKV1uFXR62voXPKTnEC1WCxPaVWApjjWgOx2Mrizn/d98iupp5XgDdWiajQs7ruUSz1zagm46imaBblKXp7rugs4Kkn+4i7HFuVyXexBM1Sa77/hqQtlY9Tq01KrzzHRfvfwyRizG+v/bQdyRhyvawfV3L+HKf5pK1eR8jKTJX36yjWB3DLvLwryf/SsfvfMnuKwuphRM4YsXfBGAifPVPlX5ZW7177n/CEBBYTjTffVX40LeOPwaM+vTM7scqhox7oKSU77g5pe78RW7MLHQlT8Fm1+HwBF2NawlN67e9+5cC3ZXf0UzUlWAp1dd3FtffILVTa8Q0nVmBXyZGVfe4oHrt+Tkq6ATs+eR257AE8wDwElkwP08DhX4cuKFjIqlq6vxXnKnDhzTd7SLPrOA6gY1g7F2r+pyHmvsYdaysVhsOvPeO1q1tXQepvMSAEbXr8LUoGrBpQOfP9/J3OU1TLuogmmXDB6uCitzKa72YKRMXlypKpczLqlk7OxT3aFOCDHcRmTQSXWrKkc8vZicNRVlY8kk5lTnUXIGBiKDGpR96Q2TuXnFhSy8bixjZxfz3s/PpGDaaIrTF6/JbYtoi/WAt4qG8DRiYTWmI5RTQfN71d5feZ2dPJubQ14ILElI2NSn3qoy9feChXNxRPorOl16CucPv0fdhz6Md8tB8nrSS+y7WtEscWZF47zvcIwep1roMeIuYXQrHAkc5JqpBZnZPZ3eKpI2VRkqHVXJ1dPLwKM+8YbdpeSHTDRLlLbwIQiqdheMH039ZhVOpl9cxfhvfwXPsmXUXnETq1Oz+HXqaqaNKqJoouoOiNl9jNdVd8Eeuw2/oQagBpzt3HHBXTTmqTE64ZiFRDy9nlD3IZqtVkrCaqxFmWfgHinTFyyn5pAKdBOaFtC7dRf11WoDxgkXFvHipEfATJG0ueletwYeu4Gvd/0gM+Mqr+j4hSdL0xWy1jr1Wrpmz8JaWooRClH39Dr2bFcX3FnaZmwOK5qmseSDqgLV3aK+N2F+GbnTJjG1cCqrPryK3773tzjSaxSV1Hj5yNfmcf3dF6jQUjoVqhdTWNQ/TsfRXsbMOhuOBMSqJlO7T40vm7ig7Lj2noimaZnuq46iGRR0a0Q0jYawlaQ9PT6nfGA3lD5uLL70WKymrQf5c676P3OJMYGIU13s+2am9emr6MQceYxqh4KQev1ctoEb7xR4VNBxMJbpTWpmkyPegVY+9YTnkLtoIbMqOhl/4I9oZori9i3MX9r/Myup8VI9OQ80nZQ1B2e0k+L2LTQUw8ya42dxLrpuHJfeMPmkC99NXtzfRVVS4+HCD0844X2FEOeeER10Ik71SXTjkst4dtxFLJt66heNU+X22pm3fDTLPzeDvFI37tmzKW9+FYCJ7fMIRSJEx1zE/qjaTsEdVuNUDvZWEHEUEPVbSRka8+OVJGzp7hVMqivUBWn2qBoCjv6KzphQHYEnngDTpCsXHHFVhbClXGimxrc6OokEqwjmqqATdeQzrtnkoAW+Xb2ZBTkqqLTiycxC8lYV4LRZuO8f1ZTcsLuUwpCqLBTm7M1MLc+pmUtLbQBNg/HzSnBfcAFVD/yIZe+Zyy2Jr/JD7SbGFefgKsnHke7CuyA+E4A9NgcBi3r9w75erhp7NT2lEawJdTH0t6UXXOw5RIeRj8McBaZBZc3AYFo0dTZFnTvICTVCXGfNc9305E9Ew2Dx8knMLJuBnkyfY2gM7PkzHSGIp4NOcdHx3UClY9SFtLMpRDyaRNN1vFdfjaFZeXW1al9583qqJuZlHlNc7WHiwv7q4OTF/e8tn8OXCTl9Smq8/UsSAMz/NDklMXLCzTgj7VhNG1ccUYOqa6ffgJEyqZ5WmNkn61Rlgk7hdEq6deosVsJdtqPG5wwciO+ZNJWCblXJqItfwKuuHDBNJrzVTdSlKnLHBp2cfPUziTnyGdVh4o2q18+dM/DXTWG5ncrGtQBUhdSAfWeqHQpOPoup5O4vUX3kZS5a92Wm7/wl3ssuGfD9+e/vf/yowy+hmwZ7p3gpcZ94Cv7JTJxfisNtxZFj5arPTMdiG5G/NoU4b43I/7GpLjVmImZXv9Q3+9WnuTM1PudkXLNmUdC9F2e8G0fKzbTWC2kqnUNdTK0oO3X3b8h3B0ildPZMvA6LoVHdBh88VJLptorq4HOri6LdaqclP5UJSLPr1RTr3js+zu8u07HHVBXCG/fyoQ4n4xMJ/P7yzFTimCOPMS0atXYbeW/cT3XqMADtITtmei2Y3Dx1IStJz0iJ2704oy4wTS6Z1k1uemp5zKK6HCon5ZPj67+QXzKxmLuvnMh9H56J1aKjaRq56YGoVUEVuPZUvw9Ts6EZSSbOmIVNt2FU5uAOq1DS3agGM9N9CCOspmp7gofxThjY5WCrqkJ32Bld/xwAzQH1Mx5THMJT4OTy6suJ6mrwdqdtHoy9lKayKzOhzpM/8KINkONz4Clwgglt9So4eq95L4eql9FrerAbYcYffBLn1CkDHrfounG4PDbKx/syM55O2fQPYbvhR1gL8phwUC0E6XcuYu+Ej9EcLUDTNZZ+ZPzQjgmUjfXhcFtJ2nKJ5IziYMxJTrvluBlXfUqnzcPnP4ieihA3fOSHKvlE92SiB9swdBu6Drn5A0Obp0AFnairkJIeC7kx9drm5A28n61mDOMPPoUr2p65rUDvAvvxWxkczTl1Kt5rrsGaiuEYMwb7MXvclY31MfOyKkZNyeeFOyzceauFtk9cNoRX6Zjny7XxD/+6kE/82yK8Rce/P4QQ57YRGXSS6YpOwqZ+ofpNjXHFOYwrzj3Zw84Ix9Sp6DYrNekL8cKGa9myexxJ04kr2oan9xC5Fa8A0JOnFk37sH8Bxa/szgwWTdgGltk7y3KYuvs3VNX+kYqmTVjy8tixsIS9VVqmolMRyeUbwQOYBrS3OUFL/+g1nepOL7UuD4Q7IBEC3UZ3t+pmsCRDmU+wdqeVnPT6QFFnGa4YrOveTEmP2gSxpUW9fsd2p+i6xhcvn8D7Z/WvjeKxqC6nnF41kHVbs7rYuaIdLF2g1lDKLy/DEVPdVx17VNUq1XMIV0h1HeT17Mc+duyA59IsFnIvvoiS9i3kJFSYwjSYfYm6kL9v7PvozlGhsD1VDjc9Q+/22kzXVd/5HevYcTrhorHUj1az6CbseQxbMoxjysCg4ylwctN/LOH6L815201tj6NpaBfcSM6ixRR3bGNM62oAGivVgngzLqk8rvpyKnSLTvn4PAB6fOPZEXMxusUk7FLB97igM3UuumlQ2KUGAY/uns4HVkeJpLd+yC10HTdGyFvkxJFjVVs25FTiiavnyz3m/5d78VLc3hBTdv0ms/OxN+fUdjEu+cqXyb3sMoq/dNeg37/oYxN5/51z+Pb7/pPPX/ttvjT//53ScU8kJ88xsOImhDhvjMigk+hUFZ2URV0oklb41NKxJ3vIGaPb7TinTqWiaR1J134sppWWDerHUNqyCU2DI/pLAJgWD0mLg5mrdmCEQoSK0uvjOAcOmI1WF+HtPcTEhpfRgLyPfIQ94TrafaA51WDkqOlDN00iTCZgGziQ0hcpoBE7mWX1CscT7lBdRTYzOOC+eemLq5piDlpvDGcCenOr6OlMYLHqjJ3z9gM1fS41+NTsTU9hT6/J4g634Rmrgszo/HGErekNHetVaOkMNFAWUF0TeT0HsR+z3xpA2b33YisrZcw+VQkpb1lPySJVbSp0FeIZr7qngikPuxs2U763k7ijL+g4jjseDBynY6QM1jyyF1OzUti5g5K2TWh2O45B2mK1W97R7BzXPLU20ujdf6SwU81Kc+RYmf++45/rVJWPV+fSkzeOUMDO6FYrkWMWC+xj9Xjw+2wUdqktKBZ0LcDce5CoTw0wP7bbCtRYoLL06+X3jiFhVSHRUzGwW1CfeBmjP5ZPdf4OJh74I76eA4yuDB53vMHYSksZ9dBP8aYHhp+Iy+riAxM+QLFbBg8LMVKNyKDT06I2ETQs6pf6n+6+hE8srD7ZQ84o1+xZaEBl5DkCjs7M7aVtG7F5TOodfiJW9Qs/6izEFlLdNq1j1fon1pyB6/yYo49az8NiIf8fPs7+7v2gaVinVIFpYKITMXyE4xMz43P6xB35FDZH6aiaq24onkSsOz0LRh84UyY//Yk/7C5ldG+K0h51e+MotfPz6JmFOFwD2zeY/Hz11gsE7RTYCyjzqwuRO9GNpUBVeWqKptGZrr50tSchHqIhYpAfVRWjAtqx5h8/psZaUEDVD39ASc8OFr/+r8xIvIklt7+asPRSdXGMuEp56affxJ6CqHPw7pU+pX0VnboA6586SPMBP1abxqR9j6EBjokT0WxnfkPGnPnqddUwmdH6LDMvreCqT08fdCr0qapIV3T8vnHM3W8QyZ2Iodtw++x4C48PLsHKPAq61OrMqVAecVsOqbmqK8h3gq6cvnFNAe+Y/m7BmmPGwNlc6J99noqPTGJ+5dO8p+c/KLps3rGHEkKId2REBp3elnYMzYKmqbEErtyzu2OwK72eR9WhDl6csBLTmqLY2UVOuAWnN0q9zUZvOgBFXCoAaDYbTXkqjLk9Ay/G1qpKEukij+fyy6G0mHp/PQB5cxdgTw/oDRv5hI4k6c1VlaG+acRRpxqQXH/h7TDhKlh8O6mQ6mqxp5fz75NX1h90pvUkKOlRdaDuAjWoeNwFpzbgs6jMgTURJp6yskC7lIKwepzHncp084ypWEBt0UE0I0Ug6aZ7fy2HkmpGTk6omdzRFSc6PK7Zsym952u4op14L1464HsLp8xDM2IYuo1Jr4YwNJ1EenkB9wm6ropHebBYdaLBBFtfVOOY3nPzVLwVeQA4j+m2OlPs48djSYe54g++j4s+PplRU97ZppDF1R4suknClsuo7lLai9T7ccys4szGoEcrmLMAZ6yHnGAjaDp1kz7MgYD6eRXXeI67PzCgotM30NtbPsjAaXcB2s3Pknf9+6m43I5l+vJ3dG5CCHGsERl0Qu1tmTV0wDylCsSZ1Bd0vA1d9DgPUf/+F1lsVzuK231J6m02/E5VderMVV0EnmXLCMfUj8t3zODP/Nwi9lZqmJpGwS23UOuvJWkm8dg8lC68ODMguTMxi+Dm7YRyVUDo2x8p5ihgbItJnZ6CG34PoxZgRtUF3+0eeOE7uqIzPpCipAcizkKithI0XaN66qldhO0lRZmumCuMD1BqjlbnVtAfNEqKp7K3Okx+t6om7F1TR0dMBQqf/wC5F1100uco+MQnGL/6ZUr+eeD4DN2i49LVa2KlbzabhqZruDyDBx2LVae4ur8qNO+9o5kwr5TCT38KbDa8y68+pfMeKk3TKPr853DPm0fBJz95Ro5pseqUVKjz7MmbQEeRCqljZxcNev9ZX7qX8hUrMu+XxqIFGCmTcReUDJh6fbSSdEUn6irCsKTfSycIkdic8MFfwD/vg4LT75ITQojBjMigk+juzAxEdjgZ9FPsu8laUYG1uBgtZTCmBVqSjSQPHgQgVpAiomsEnWocUX3hBOIWGwW33IIRUQOEC4/pXihwFvCDD+j88RtLcF8wh/09agfwCfkTcE2biiOpur6at7sJu0owdBs2p4XK9HToqDOfsc0mdX61XoppmpBUr4/HOzBUZSo6rmKqeooo8Zt0FqrxL+XjfDjcp1YdsxYVUdyhNrRs3Rkmkd5XK7+qf3aSrlvw5dpIpNS2Dvt3poiH1QDtvJ6DmYX7TsZWVoZmPT7IFheq8wrllBPzpbvNPLaTrqdSOUlVVkbPLMpsF5D3oQ8x+a1t5CxZ8rZtOV0FN91Eze9+i7Vo8CByOiqmqIrM4arLiNu92KwmlRMHX2HZkptD3geuZ+IHF2duq5qcz5X/OPWEr5fDZSW/pH/qv40YVtvb7As11AHbQghxCkZk0NH9PSSsalCt0312qzmgPqW75qjxNrNrDdqDrcQOqE09Wy9VM44sPjULpTF/LB9d/i1SEydjTajbykoHzrbJd+bT69Y4WKjG1ezv7g86mt2OKz2mJ2I4iFSralJRVS6ewvQ0YEcB5V3Q2Kba0BxqRiO91URx3oDn8uQ7sVjA1K0ka10s3GfSmd6Z++120D6ataiIgu7d6EaC3s4oJhp6KoZ39MAKwWibl21VO9CNBIGUDz2mBo07c7pxjD39T/8lo1XlKZhTgWWWuoCfaHxOn7lXj+aa22Zy9WemDwjHQ55RdQ6onKrCXditxs1UT1RdcydTNs5H1eR8Rk0tYPnnZrztejJlE/qD07GLBQohxNkyIoOOI9RLIr0qsiNneKaM9u2XtHi3id7YhhmLoTmd1M5R3TGeIhVC8k0LMauDTYe6cRvqglpeMnBmTKFTBYyu9G7Q+7r3ATAhT81e8pSm199xFmJc/iEAikZ5jlrvpAAdiKSrSru7doOeHldRMnD9F03XyEuvyhx2l+KK2ehOP0/NjFMPOo4pU7A7reR37cnc5o6046gZuCZKTU4566fEKEh3cwE4ol3YLpp8ys81mJLp6nl6PdV0lKqFEN0nmHHVx+awMHpmUVYsGFc6xpuZ0g0wbknN2z7GYtG57q45vP+Ls7E73/4DQtnY/jE5x3aBCiHE2XL+/8YeIsMwyIlGSKaDjtNzZrZ8GCrPey5Dczqp6IaL3lKVGMf48RwKqYGupeXq07A3paGZcMuv3sSVnv/tyR/Y5gKnqk70BZ2jKzoABReocS3N5UvYtUs9V/GoXHLTx0lZXCQtTqxNbcRTcfZ07iFlVRep3NLjB5D2TUHuHb+Y7vxJGBY7uXk2CspPfV0X3eEg9/L3ZLqvQE0tt1UPnP022jeOXreGbmzM3JbnP3BK3VYnUzxVjX2KOQs42KjCbl7JyFkMzu60km9XaxnpZmpI1bhT1TfzCsA7+t1fjFMIIQYz4oJOe2sXVsMgkR6MPFwVHT0nh9zLLgVg+UaVYByTJnIooHYPry5TG4LqgBcNlwk6GibHzxLrCzqBeIDOSCetYbXI3vh8tXLu5PfOYOKCUixWHdNQz1VS48XmsGSmKUed+ZR2mTQEGtjfXAuaGk+RW3n8BXBUesBxnWs2+yaqXbtrZhQPuQvHe/Vyijq3ZyoLrmg7trKBU5BnTLoegPU1u7CkYurGxEHKZxy/b9FQ5PgclFTnoGlqXZklHxz/jtamOR9VTFQ/x7Ji45QqNENVUJaDPb3mk6f81HZYF0KIM+3sD1AZZo31TbiAqEMFnXeyHsk75X3ve+n923M4VZEF58SJHAr8EYDReTUkCg38bRH+dMtCTIeFP39/M65c23EL0PkcPnRNxzAN/nfX/wJQnlOON72SsjPHxpX/NI2LPpbg4OY2rDadwvR2DrkFDqKhBFFHAWXdLdQF6uioO8IE1E7SjtLjF1qbsqScYHeMDX+uI+pQF8vRM4c+UDZn6YU4nRp5/gP05E0kzxE5buDwmMoFTMqfxGuT97Dsyb/QUTiL2rKdeB1D3FJhEB/62gKS8dS7cpE/H8y/aQGppw8y+8p3Zw0pTdcoG+ujYVcXuQXDUzkVQogRV9FpOZzeEdvZN0Zn+C5yuRdfTOyoVY4t48dypFftRD7aOzqzGFuwK4YZVOvZDDb9Wdd08hx5ADy842EAphdNP+5+zhwb0y6qZNKi/gG/mXE6zgLKuk22tG0htzHdlWZG0F2Dr3y74H1juPjjE0FT6/H0zUgaCt1ux3P55UzZ/Vsm7/kdowojg95v+ZjlhFwatUVrmbLrfvYuPjPVAV3XRmzIAVXVuvzmqRRWvHtbnyz+4DguuKqaCfOl60oIMTxGXNDpalRbCmQqOqc4HfrdoDscNM/t/zTdVekhaSZxWV2UuEvwppfXD7RHOLhZbW5ZMX7w3aqnFaqZT6M8o/jC7C/wb4v/7ZTa0PdJO+bIp7wLnq9/nvIe9Rwuh3myhzLj0io+8rV5fPDLF2Czv83U4RPwvnc5rlgXFS3rcVSPGvQ+y8eoReQeWmZwy5cs2GvO3irW4p0pqvKw+APjz/paVUII0WfE/fbpbVVBJ943GHkYKzoAvZfMhlfrCOW7aLOoRexqvDVompbZKbmzKcSRvWoj0kkLywY9zvcv/j5t4TbG+MYMaaxM38DmqLMAbwSCXa0UhKaB6yQLvB2lpOaddSHlLFqE7vNh+P3HDUTuU5FbwZySOWxp2wKobjkhhBDiVIy4ik6kQ22tkLL2dV0NX0UHwLZoHj9frvOXmyZSH6gHVNCB/n2EDm3vIBlL4S1yUjZu8IpOrj2XsXljhzwgOLdATamO5qqxOKXd4I2q6eg5RYMv738maXY7hZ/+FJaCAnIvueSE9+ur6oAEHSGEEKduxAWdZLeqjKQswz8YGaA0p4y/z9b5k+cAK3euBNT4HABvsaq2mOkepIkLy8744nR9iwbGXGp2VUWXiSOZ3m161Jlbifdkij7zGSa+9uqgu3/3WVazDEt6JlhZzuBVLSGEEOJYIyroBGNJ7L0B9YWWnl4+DCsjH21y4WQ8Ng+xVIyOiNrfanKBWgzPe8zO0JMWnPkLfN9g5Jieg6HpTD9kZtbQ8dacOwNIC12FXDvuWlxWF3NL5w53c4QQQpwnRtQYnYbOMN54CBMts3P5qe7N9G4pcBaw6iOrONhzkCO9RzAxeU/1ewC1qJvLYyPSm6B0jJe8UvfbHG3o3B47ukXDSOnE7T7m7e9m97T0Tt4FZ/753olvLf4W/7r4X7Hpw/szE0IIcf4YWUGnK4Q3HiJpdaKli1nnwmyQHFsOM4tnMrN45nHfyyt1E+n1n3AQ8jul6WrQc09rmLCrmIKebuLp9XdyTmEw8tlk0S1YOL3ZXUIIIUamEdV11dgTxRcLkrSqLiGLlXN+36KLPjqRxR8Yx9SLKt615/Cltz6IuEpI6TaS6Z3d3d5zK+gIIYQQQzX85YyzqCMYY0Y8RNKhBt46XOd2yAEorvZQXP3uzn7KK3ZziE7C7mLidvVcuq4WAhRCCCHOZ+f+lf4M6ugJ40lEMhUdh1sqFnBURcdXlem2cnusZ3yGlxBCCHG2jaigE+pQu3vHh3lDz3NNXol6PSI5Zf3jc86xgchCCCHE6RhRQSeaXiwwlN6/ye6S2TtwVEXH4iXmUPtIyfgcIYQQ2WBEBZ1kl6rohF19XVcyBgXUfle6VcNAJ+BRqzLn+BzD3CohhBDinRsxQccwTEy/2ksq7Dw3Fgs8V+i6ltluonf0PODU9rkSQgghznUjJuh0h+PkxkIAxNPTp8+FNXTOFb70OJ1QVK1TI11XQgghssGICTodwTj2VAIgM+vKLhWdjL5xOn3c0nUlhBAiC4yYoNPeG8sEHaNvQ89h3v7hXNI386qPVHSEEEJkgxETdDqCRwUdvW/WlVR0+hxX0ZGgI4QQIguMmKDT3hvDbiTVF7oMRj6WVHSEEEJkoxETdI6u6KDJ9PJj5eY5Mvt+OXNsWKwj5q0hhBAii42Yq9mAig7SdXUsTdfwFavXRaaWCyGEyBYjJ+ikKzqGpqNrTkAGIx8rE3Sk20oIIUSWGDlBpzeGzUhkppYD2F2WYWzRuadvnI5UdIQQQmSLERN0OoJxclIRkukNPa0OHd0yYk7/lExcWEbZWB+TF5cPd1OEEEKIM2JEDFJJGSZdoRguM5Kp6Dik2+o4RVW5fOgrc4e7GUIIIcQZMyJKGl2hOIYJTiOWCTpOmXElhBBCZL0REXTae2MAuMx4putKZlwJIYQQ2W9EBJ2OoAo6DiNBQrquhBBCiBFjRASdvoqONZXKVHRk53IhhBAi+42IoNNX0bGmDNm5XAghhBhBRkTQ6a/oGEfNupKgI4QQQmS7ERF0+io6etKQrishhBBiBBkRQac9HXQsSbM/6EhFRwghhMh6IyLodPTGsRgpdJP+WVcumXUlhBBCZLsREXTaeqPYUwkAGaMjhBBCjCBZH3SiiRTd4QR2IwnQv2CgBB0hhBAi62V90GkLqPE5uVoKOKqiI4ORhRBCiKyX9UGnJRAFoNylk9KtGBY7IF1XQgghxEiQ9UGntS/oOHVSFlf6VhO7U4KOEEIIke1GTNApcWokbGp8ju4w0XRtOJslhBBCiLNgyEFn7dq1XHvttVRUVKBpGk8//fSA77e2tnLLLbdQUVGB2+3m6quvZv/+/Sc95sqVK9E07bg/0Wh0qM07Tos/HXTsqcz4HIvjHR9WCCGEEOeBIQedUCjErFmzePDBB4/7nmmaXH/99dTW1vLMM8+wZcsWampquOKKKwiFQic9rtfrpbm5ecAfp9M51OYdp2+MTqk1nplxZXVlfSFLCCGEEMCQB6osX76c5cuXD/q9/fv38/rrr7Njxw6mTZsGwE9/+lNKSkp49NFH+fSnP33C42qaRllZ2VCb87b6uq4KtWimomOTGVdCCCHEiHBGSxuxmJrKfXQlxmKxYLfbWbdu3UkfGwwGqampoaqqive9731s2bLlbZ8rEAgM+DOY1vT08nyODjqWUz4nIYQQQpy/zmjQmTx5MjU1Ndxzzz10d3cTj8f53ve+R0tLC83NzSd93MqVK3n22Wd59NFHcTqdXHjhhScd27NixQp8Pl/mz6hRo467j2mama4rrxElpaup5XaHVHSEEEKIkeCMBh2bzcYTTzzBvn37KCgowO12s3r1apYvX47FcuIqyqJFi/jkJz/JrFmzuOiii/j973/PxIkT+fGPf3zCx9xzzz34/f7Mn8OHDx93n55wgnjSAMBlREhZ+oKOVHSEEEKIkeCMlzbmzp3L1q1b8fv9xONxiouLWbhwIfPmzTvlY+i6zvz5809a0XE4HDgcJ58+1VfNKcixo0dCGJmKjmzoKYQQQowE79r0I5/PR3FxMfv372fjxo1cd911p/xY0zTZunUr5eXl76gNfUGnxOMgGQ6SSs8rdzgl6AghhBAjwZArOsFgkAMHDmS+rqurY+vWrRQUFFBdXc0f/vAHiouLqa6uZvv27dx5551cf/31LFu2LPOYm266icrKSlasWAHAvffey6JFi5gwYQKBQIAHHniArVu38pOf/OQdnVxbOuiU+Zwkm4OkLCo4OR32d3RcIYQQQpwfhhx0Nm7cyGWXXZb5+u677wbg5ptvZuXKlTQ3N3P33XfT2tpKeXk5N910E9/85jcHHKOhoQFd7y8m9fT0cOutt9LS0oLP52POnDmsXbuWBQsWnO55AdDiVzOuyrxOEpEwKT1d0ZGgI4QQQowImmma5nA34kwIBAL4fD78fj9erxeAe57czqNvNnDn5RN4/x++xMauK+komsmlN0xi2kWVw9xiIYQQQgx2/T6TsnqJ4Naju65iscysK6tdZl0JIYQQI8GICDqlXgdGPJZZR8cmQUcIIYQYEUZI0HGSisUx+io6jqw+bSGEEEKkZe0VP5406AjGATUY2YgnMxUd6boSQgghRoasDTptvaqaY7foFOTYMRPJzBgd6boSQgghRoasDTp93VYlXgeapmEmUv1dV/asPW0hhBBCHCVrr/h9u5aXetVO6mbS6B+MLHtdCSGEECNC1gadjqAKOsW5apFAIwmmrtZHlDE6QgghxMiQtUGnMz0QuSBXVXFMo381ZBmjI4QQQowMWRt0ukIq6BTm2MEwMMx04MFAt2rD2TQhhBBCnCVZH3QKcuyQCGGiurCwpNA0CTpCCCHESJC1QaczpMboFOY6IBbEJN11ZTWGsVVCCCGEOJuyN+gE+7uuzFgQM911pVmzYg9TIYQQQpyCrA06A7quwj0YFtV1pdkk6AghhBAjRVYGHcMw6Q73V3SM3u7MGjrYZHyOEEIIMVJkZdDpiSQw0oWb/Bw7ZsifWRXZ4pCgI4QQQowUWRl0utIDkb1OKzaLjhnyZyo6ulR0hBBCiBEjK4NOZiBy36rIwZ7Mhp6yz5UQQggxcmTlVX/AQGTADPdmgo5Fgo4QQggxYmTlVb/zuKATxOjb0FO2fxBCCCFGjKwMOn0VnaL0PldGuJdUenq57FwuhBBCjBxZHXQyFZ1If0XH7rAOW7uEEEIIcXZlZdDpCKpZVwU56cHIkUhmjI7dYRu2dgkhhBDi7MrKoDNg53LAjIYzQcchQUcIIYQYMbI66GS6rqKRzDo6EnSEEEKIkSMrg86xs66MaCSzMrLT6Ri2dgkhhBDi7Mq6oGOaJt19XVfpWVfJWCxT0XG67MPWNiGEEEKcXVkXdAKRJMn0Rld9FZ1EPJ4Zo+OSio4QQggxYmRd0OkKqxlXuQ4rDqtaMycZT2S6rlxO57C1TQghhBBnV9YFnWO7rTBNkolkputKVkYWQgghRo6sCzpd4QTQ321FIkwyRabrSlZGFkIIIUaOrAs63cesoUM8RCqlY6S3gLBKRUcIIYQYMbIv6IQHTi0n1kvC7J9pZZXdy4UQQogRI+uu+v1r6KRnV8WDJIz+mVZS0RFCCCFGjqwLOoN1XSVNFXRMEui6NlxNE0IIIcRZln1B59jByLEgKTM9pVxPDlOrhBBCCDEcsi7oZPa56pteHu8l1VfRkaAjhBBCjChZF3Ra/REAyrzpKk48hNE3GFlPDVOrhBBCCDEcsi7o9ERU1WZUgVvdEAuSIj0Y2SJBRwghhBhJsi7ogBqfk+uwqi/iQcy+oGM1hq9RQgghhDjrsjLoZKo5APEgBqrrSrOaw9QiIYQQQgyH7Aw6+a7+L2JBTE1VdDS7BB0hhBBiJMnKoFN9bEVHS1d0bLKGjhBCCDGSZH3QMaO9mOmuK90hQUcIIYQYSbIy6AwYoxMNkkpv6KnLzuVCCCHEiJKVQefoio4R6cXQVUXH4pSgI4QQQowkWRd0LLpGuc+Z+doMB0lZVNCxOq3D1SwhhBBCDIOsCzplPgdWS/9pmZEwKV2CjhBCCDESZV3QqcpzD/jajIYw0hUdu12CjhBCCDGSZF/QOXoNHdPEiIQzXVcOp22YWiWEEEKI4ZB1Qafy6KCTjGGmzEzXld0hFR0hhBBiJMm6oFOVP3CxQDOlZbqunA7HMLVKCCGEEMMhC4PO0ds/9GKm6K/oSNeVEEIIMaJkX9DJOyroxEMYKS0zRsclFR0hhBBiRMm6oJOfY+//Ih7ENDRMXY3NcTmdJ3iUEEIIIbJR1gUdTTtqP6tYECOlYeiqy0oqOkIIIcTIknVBZ4B4kLihYaQrOm6n620eIIQQQohskvVBJ2ZYMTW1x5XLIV1XQgghxEiS3UEnFiRm9I/ZcTjsJ7mzEEIIIbJNdgedeJD4UUHHYs3u0xVCCCHEQEO+8q9du5Zrr72WiooKNE3j6aefHvD91tZWbrnlFioqKnC73Vx99dXs37//bY/7xBNPMHXqVBwOB1OnTuWpp54aatOOFw8SM9NBx0yi69rJ7y+EEEKIrDLkoBMKhZg1axYPPvjgcd8zTZPrr7+e2tpannnmGbZs2UJNTQ1XXHEFoVDohMdcv349H/vYx7jxxhvZtm0bN954Ix/96Ed54403htq8gWJBkkbfIoHJd3YsIYQQQpx3NNM0zdN+sKbx1FNPcf311wOwb98+Jk2axI4dO5g2bRoAqVSKkpIS7rvvPj796U8PepyPfexjBAIB/va3v2Vuu/rqq8nPz+fRRx89pbYEAgF8Ph9+vx+v16tufOrzvPZ/29hS9B0gyO0/e//pnqoQQggh3gWDXr/PoDM6aCUWiwHgPGphPovFgt1uZ926dSd83Pr161m2bNmA26666ipee+21kz5XIBAY8Oc48V6S6a4rk9RQTkUIIYQQWeCMBp3JkydTU1PDPffcQ3d3N/F4nO9973u0tLTQ3Nx8wse1tLRQWlo64LbS0lJaWlpO+JgVK1bg8/kyf0aNGnX8nWJBUunByJouQUcIIYQYac5o0LHZbDzxxBPs27ePgoIC3G43q1evZvny5VgslpM+dsCKxqjxPsfedrR77rkHv9+f+XP48OHj7xQPkSI9Rkczhnw+QgghhDi/Wc/0AefOncvWrVvx+/3E43GKi4tZuHAh8+bNO+FjysrKjqvetLW1HVflOZrD4cDxdls6xIOkzDL1b6noCCGEECPOu7awjM/no7i4mP3797Nx40auu+66E9538eLFrFq1asBtL7zwAkuWLHlnjYgFMcx0RUc/7THXQgghhDhPDbmiEwwGOXDgQObruro6tm7dSkFBAdXV1fzhD3+guLiY6upqtm/fzp133sn1118/YLDxTTfdRGVlJStWrADgzjvv5OKLL+a+++7juuuu45lnnuHFF1886QDmUxLvxUgPRtYsEnSEEEKIkWbIQWfjxo1cdtllma/vvvtuAG6++WZWrlxJc3Mzd999N62trZSXl3PTTTfxzW9+c8AxGhoa0PX+YtKSJUt47LHH+Jd/+Re++c1vMm7cOB5//HEWLlx4uuelxEOY6VPUrBJ0hBBCiJHmHa2jcy45bh5+Mg7fKeZPmz5KQ+U/YMlv4nMrPjnczRRCCCHEUc6rdXTOKfEgAGZ61pVuk+0fhBBCiJFmBAWd7D1VIYQQQgwue6/+MQk6QgghxEiXvVf/dEWHdNCx2LP3VIUQQggxuOy9+h/TdWVxnPG1EYUQQghxjsveoBMLYpqApgKO1Wkb3vYIIYQQ4qzL4qDTi2mAoaugY3Pah7lBQgghhDjbsjfoRLowUxqGrio5Nufb7IslhBBCiKyTvUEn3Ilp9AcduwQdIYQQYsTJ4qDTV9FRXVcO6boSQgghRpwsDjqdGCkyFR2HXYKOEEIIMdJkb9CJdA/ounJK0BFCCCFGnOwNOuFOkkd1XTlljI4QQggx4mRx0Okibuj9QccuQUcIIYQYabIz6JgmRLqIHdV15XI6h7lRQgghhDjbsjPoxAJgJEkYev86OjbZAkIIIYQYabIz6IQ7AYhrjkzXlVV2LxdCCCFGnOy8+oe7AUjo7kxFxyJBRwghhBhxsvPqn6noOPuDjjU7T1UIIYQQJ5adV/9IFwAJLSdzk3RdCSGEECNPdl790xWdBK7MTdJ1JYQQQow82Xn1D6crOmb/lHLdog1Xa4QQQggxTLI06KiKTspQ43MwE2iaBB0hhBBipMnOoJMeo5M01f5WmpYaztYIIYQQYphkZ9BJd12ljPQigRJ0hBBCiBEpy4OOBZCKjhBCCDFSZWnQUWN0jFRf0DGGszVCCCGEGCbZF3TSG3oCmH1BR5egI4QQQoxE2Rd04kFIxQEwDTXTSrOYw9kiIYQQQgyT7As66fE5WJ2YSXV6mmUY2yOEEEKIYWMd7gaccRG1oSeuAkxDBR09HXQMwyAejw9Tw959NpsNi0VSnRBCCNEne4OOuxBSqutKt0E8Hqeurg7DyO7xOnl5eZSVlckCiUIIIQRZHXTyIV3RsVh1mpubsVgsjBo1Cl3Pvh470zQJh8O0tbUBUF5ePswtEkIIIYZfFgad9BgddyGmqbpx7Lk2wuEwFRUVuN3uYWzcu8vlUpuYtrW1UVJSIt1YQgghRrzsK22E+8fokA46NqfKc3a7fbhaddb0BblEIjHMLRFCCCGGX/YFnaMqOpmg41BBZySMWxkJ5yiEEEKcqiwMOn1jdAro65mz2LKvh04IIYQQby8Lg85RFZ100LE6JOgIIYQQI1H2BZ30GB3TngeaDQCb4/wem/PTn/6UMWPG4HQ6mTt3Lq+88spwN0kIIYQ4L2Rf0OmqBcDMrcDQ00HH6RjOFr0jjz/+OHfddRff+MY32LJlCxdddBHLly+noaFhuJsmhBBCnPOyL+ikYuDwknCVYeiqy8rhcA5zo07f/fffz6c+9Sk+/elPM2XKFH74wx8yatQoHnrooeFumhBCCHHOy87BK+WziMVCmYqOw2kH+rd+ME2TSCI1LE1z2SynPDMqHo+zadMmvva1rw24fdmyZbz22mvvRvOEEEKIrJK1QScS9mcqOmr9nP6gE0mkmPqvzw9L03b9+1W47af2snd0dJBKpSgtLR1we2lpKS0tLe9G84QQQoiskn1dVwAVc4hFgpmKjsV6fp/msRUg0zRlvRwhhBDiFGRnRadiDvEDtZmKjsWqc3RHlctmYde/XzUsTXPZTn1bhqKiIiwWy3HVm7a2tuOqPEIIIYQ4XvYFHbsX8scQj24/qqIzsPqhadopdx8NJ7vdzty5c1m1ahUf+MAHMrevWrWK6667bhhbJoQQQpwfzv2r/VCVTQNdJxHtH4ysWy0wPGOP37G7776bG2+8kXnz5rF48WJ+8Ytf0NDQwOc+97nhbpoQQghxzsvCoDMTID1GpwgAq1U/b4POxz72MTo7O/n3f/93mpubmT59On/961+pqakZ7qYJIYQQ57zsCzrlKugkI+GjKjoaxIazUe/Mbbfdxm233TbczRBCCCHOO+f3dKTBlM0AIBmNYGj9g5GFEEIIMfJkXwLIHwNAPIumlwshhBDi9GRfAkivLxMJ9GDqaiq3xSJrzgghhBAjUfYFnbRYbyjzb92WtacphBBCiJPI2gQQD4Uz/5aKjhBCCDEyZW3QSYX79rYy0C1Ze5pCCCGEOImsTQBGJAGArhnD3BIhhBBCDJesDTpE1QqBllPfWkoIIYQQWSYrg45hGuhxE5CgI4QQQoxkWRl0/DE/tpRaLNAqa+gIIYQQI9aQU8DatWu59tprqaioQNM0nn766QHfDwaDfOELX6CqqgqXy8WUKVN46KGHTnrMlStXomnacX+i0ehQmwdAd7QbW0otFmi1nt8lnbd7vYUQQghxYkMOOqFQiFmzZvHggw8O+v0vfelLPPfcc/zud79j9+7dfOlLX+KOO+7gmWeeOelxvV4vzc3NA/44nc6hNg+AzmhnJuhYzvM1dN7u9RZCCCHEiQ15U8/ly5ezfPnyE35//fr13HzzzVx66aUA3Hrrrfz85z9n48aNXHfddSd8nKZplJWVDbU5g+qKdmE10l1X9vO7ovN2r7cQQgghTuyMlzuWLl3Ks88+S2NjI6Zp8vLLL7Nv3z6uuuqqkz4uGAxSU1NDVVUV73vf+9iyZctJ7x+LxQgEAgP+9OmKdmHtq+gMFnRME+Kh4fljmkN/UYUQQghxWoZc0Xk7DzzwAJ/5zGeoqqrCarWi6zq//OUvWbp06QkfM3nyZFauXMmMGTMIBAL86Ec/4sILL2Tbtm1MmDBh0MesWLGCe++9d9DvdYU7KcEOgM1pO/4OiTB8t2LoJ3cmfL0J7DnD89xCCCHECPOuBJ3XX3+dZ599lpqaGtauXcttt91GeXk5V1xxxaCPWbRoEYsWLcp8feGFF3LBBRfw4x//mAceeGDQx9xzzz3cfffdma8DgQCjRo0CwO9vI2VxAGBzDxJ0hBBCCDEinNGgE4lE+PrXv85TTz3FNddcA8DMmTPZunUr//Vf/3XCoHMsXdeZP38++/fvP+F9HA4HDodj0O8F/e0kLWogs81lP/4ONreqrAwHm3t4nlcIIYQYgc5o0EkkEiQSCXR94NAfi8WCYZz6VgymabJ161ZmzJhxWu0IB7ow9NEA2B2DnKKmSfeREEIIMQIMOegEg0EOHDiQ+bquro6tW7dSUFBAdXU1l1xyCV/+8pdxuVzU1NSwZs0afvOb33D//fdnHnPTTTdRWVnJihUrALj33ntZtGgREyZMIBAI8MADD7B161Z+8pOfnNZJRXq7SVkmAWBznt+zrt7u9RZCCCHEiQ056GzcuJHLLrss83XfOJmbb76ZlStX8thjj3HPPfdwww030NXVRU1NDf/xH//B5z73ucxjGhoaBlR9enp6uPXWW2lpacHn8zFnzhzWrl3LggULTuukYr09JPvG6DjO76Dzdq+3EEIIIU5syEHn0ksvxTzJFOmysjIefvjhkx5j9erVA77+wQ9+wA9+8IOhNmVQiVSCVDiUGYx8vq+j83avtxBCCCFO7PxeNngQPbEeXHEwsqSiI4QQQojTl3VBpzvajTNO1nRdCSGEEOL0ZW3QSVn6FgyUoCOEEEKMVNkXdGLdOBP0Lxh4no/REUIIIcTpy76gE+3GGTf7g450XQkhhBAjVpYGHSToCCGEECILg06sG8fRQUfG6AghhBAjVvYFnWg3zoQNNHVqMkZHCCGEGLmyL+jEunEk+zfytErXlRBCCDFiZV3Q2du9NxN0LBYTXdeGuUVCCCGEGC5ZF3SSRhIPPgBs1vM75KxYsYL58+fj8XgoKSnh+uuvZ+/evcPdLCGEEOK8kXVBB8Cr5wFgtZ/fp7dmzRpuv/12Xn/9dVatWkUymWTZsmWEQqHhbpoQQghxXhjypp7nukn5k9BS6rRs53nQee655wZ8/fDDD1NSUsKmTZu4+OKLh6lVQgghxPkj64LOJ6d+kmTiGQBszsFPzzRNIsnI2WxWhsvqQtNOr0vN7/cDUFBQcCabJIQQQmStrAs6l1VeyrbkyYNOJBlh4SMLz2azMt74xBu4be4hP840Te6++26WLl3K9OnT34WWCSGEENkn64KOFo33Lxbosr/Nvc8fX/jCF3jrrbdYt27dcDdFCCGEOG9kXdAxIuH+ncvdtkHv47K6eOMTb5zNZg147qG64447ePbZZ1m7di1VVVXvQquEEEKI7JR1QccMhzMVHfsJFgvUNO20uo/ONtM0ueOOO3jqqadYvXo1Y8aMGe4mCSGEEOeVrAs6RiSSNRt63n777TzyyCM888wzeDweWlpaAPD5fLhcQ68MCSGEECPN+T3/ehBGKJQJOuf79g8PPfQQfr+fSy+9lPLy8syfxx9/fLibJoQQQpwXsrSikx6jc54HHdM0h7sJQgghxHkt+yo6R3VdnWiMjhBCCCFGhuwLOkcNRj7fKzpCCCGEeGeyOuic72N0hBBCCPHOZF3QMaWiI4QQQoi0rAs6Rvjo6eVZN9ZaCCGEEEOQhUFHKjpCCCGEULIv6GTRgoFCCCGEeGeyLuikQiFSutrjSoKOEEIIMbJlXdBJROKgqdOSoCOEEEKMbFkXdJLRRPpfJlZb1p2eEEIIIYYg65JAPJoEwGrV0HRtmFvzzjz00EPMnDkTr9eL1+tl8eLF/O1vfxvuZgkhhBDnjewLOv4QADb7+X9qVVVVfO9732Pjxo1s3LiR97znPVx33XXs3LlzuJsmhBBCnBeybqGZZDACgM1lG+aWvHPXXnvtgK//4z/+g4ceeojXX3+dadOmDVOrhBBCiPNH9gUdPT21/CRBxzRNzEjkbDVpAM3lQtOG3qWWSqX4wx/+QCgUYvHixe9Cy4QQQojsk3VB51R2LjcjEfZeMPdsNWmASZs3obndp3z/7du3s3jxYqLRKLm5uTz11FNMnTr1XWyhEEIIkT3O/4Esx0hZ7ED2bOg5adIktm7dyuuvv87nP/95br75Znbt2jXczRJCCCHOC1lb0TnZGjqay8WkzZvOVpOOe+6hsNvtjB8/HoB58+axYcMGfvSjH/Hzn//83WieEEIIkVVGZtDRtCF1H51LTNMkFosNdzOEEEKI80LWBR0j3XWVDasif/3rX2f58uWMGjWK3t5eHnvsMVavXs1zzz033E0TQgghzgtZF3RSfbOu7Od/0GltbeXGG2+kubkZn8/HzJkzee6557jyyiuHu2lCCCHEeSHrgo7pyQPA5jz/g86vfvWr4W6CEEIIcV7LullXZm4ekB1dV0IIIYR4Z7Iu6BhuDyBBRwghhBDZGHQcOQBYs2CMjhBCCCHemewLOta3n14uhBBCiJEh64JOUlN7XEnQEUIIIUTWBZ1gbwqA3HzHMLdECCGEEMMt64JOMm5gten4Ss7PlY+FEEIIceZkXdABKKjIQde14W6GEEIIIYZZVgadoqrc4W6CEEIIIc4BWRl0Cqs8w90EIYQQQpwDsjLoZGtFZ8WKFWiaxl133TXcTRFCCCHOC1kZdAorc4a7CWfchg0b+MUvfsHMmTOHuylCCCHEeSPrgo4n34nDbRvuZpxRwWCQG264gf/5n/8hPz9/uJsjhBBCnDeybvfygoq3r+aYpkkybpyF1hzPatfRtKHNCLv99tu55ppruOKKK/jOd77zLrVMCCGEyD5ZF3ROpdsqGTf4xZ1rzkJrjnfrjy4Z0qrNjz32GJs3b2bDhg3vYquEEEKI7JR1QaegInsGIh8+fJg777yTF154AafTOdzNEUIIIc472Rd0TqGiY7Xr3PqjS85CawZ/7lO1adMm2tramDt3bua2VCrF2rVrefDBB4nFYlgssqeXEEKI/9/evQdFVb9/AH+vtKxIuIkXlhUhhtEhXGMSy2BIjQohLzg2CuUodnHGEvNWg1YMNtVgN/5oCPMPNR2b8I8EmaEsmABl1DJAQ3IUx03UQIwSRAVWeL5/9OP8PNzRZZc9vF8zOwOfc5nn4fmcOY/HPedQTzTX6BjHevS5jk6nc4mXfj7zzDOoqKhQjb388ssIDg5GcnIymxwiIqI+DPiuq8OHD2PBggUwm83Q6XTIyclRLW9qakJSUhL8/Pzg4eGBRx55BNu3b+9zv9999x1CQkJgMBgQEhKC7OzsgYYGANBp6NUPXl5esFgsqo+npyfGjh0Li8Xi7PCIiIiGvAE3Ojdv3kRoaCgyMjK6Xb5hwwYcOnQI+/btw5kzZ7BhwwasXbsWBw8e7HGfx44dQ3x8PJYvX45Tp05h+fLlWLp0KX755ZeBhkdERESk0ImI3PPGOh2ys7OxaNEiZcxisSA+Ph4pKSnKWFhYGJ5//nl88MEH3e4nPj4ejY2N+OGHH5SxmJgYjBkzBt9++22/YmlsbITRaERDQwNGjx6tWtbc3Ayr1YrAwEDNf6l3OOVKRESur7fztz3Y/YGBkZGRyM3NxZUrVyAiKCwsxLlz5zB37twetzl27Biio6NVY3PnzsXRo0d73KalpQWNjY2qDxEREdHd7N7ofPHFFwgJCYGfnx/c3d0RExODzMxMREZG9rhNbW0tfHx8VGM+Pj6ora3tcZu0tDQYjUblM2nSJLvlQERERNowKI3O8ePHkZubi9LSUnz++ed44403UFBQ0Ot2nZ8WLCK9PkF4y5YtaGhoUD6XLl2yS/xERESkHXa9vfz27dt45513kJ2djXnz5gEAHn30UZw8eRKfffYZnn322W63M5lMXa7e1NXVdbnKczeDwQCDwWC/4ImIiEhz7HpFx2azwWazYcQI9W7d3NzQ3t7zu6XCw8ORn5+vGvvpp58QERFhz/CIiIhomBnwFZ2mpiacP39e+d1qteLkyZPw9vaGv78/Zs+ejbfffhseHh4ICAhAcXEx9u7di/T0dGWbFStWYOLEiUhLSwMArFu3DrNmzcLHH3+MuLg4HDx4EAUFBSgpKbFDiv/vPm4wcxm9NZRERETDzYAbnd9++w1PP/208vvGjRsBAImJifj666+RlZWFLVu2YNmyZfjnn38QEBCAjz76CKtXr1a2qa6uVl31iYiIQFZWFt577z2kpKQgKCgI+/fvx8yZM+8nN4Ver4dOp8O1a9cwfvz4Ab893BWICFpbW3Ht2jWMGDEC7u7uzg6JiIjI6e7rOTpDSV/34Tc1NeHy5cuav6ozatQo+Pr6stEhIiKXMNjP0dHcu6568uCDD2Ly5Mmw2WzODmXQuLm54YEHHtDkFSsiIqJ7MWwaHeC/RoAvwiQiIho+7P4cHSIiIqKhgo0OERERaRYbHSIiItIszXxHp+NuKr7ck4iIyHV0nLcH665ozTQ69fX1AMCXexIREbmg+vp6GI1Gu+9XM42Ot7c3gP8eRjgYf6ihqrGxEZMmTcKlS5cG5fkDQxXzZt7DAfNm3sNBQ0MD/P39lfO4vWmm0el40rLRaBxWE6TD6NGjmfcwwryHF+Y9vAzXvDu/J9Nu+x2UvRIRERENAWx0iIiISLM00+gYDAakpqbCYDA4OxSHYt7Mezhg3sx7OGDeg5O3Zl7qSURERNSZZq7oEBEREXXGRoeIiIg0i40OERERaRYbHSIiItIsTTQ6mZmZCAwMxMiRIxEWFoYjR444OyS7SktLw+OPPw4vLy9MmDABixYtwtmzZ1XrrFy5EjqdTvV58sknnRSxfWzdurVLTiaTSVkuIti6dSvMZjM8PDwwZ84cVFZWOjFi+3j44Ye75K3T6bBmzRoA2qn14cOHsWDBApjNZuh0OuTk5KiW96e+LS0tWLt2LcaNGwdPT08sXLgQly9fdmAWA9db3jabDcnJyZg2bRo8PT1hNpuxYsUK/PXXX6p9zJkzp8scSEhIcHAmA9NXvfszr7VWbwDdHus6nQ6ffvqpso6r1bs/5yxHHt8u3+js378f69evx7vvvovy8nI89dRTiI2NRXV1tbNDs5vi4mKsWbMGx48fR35+Pu7cuYPo6GjcvHlTtV5MTAxqamqUz/fff++kiO1n6tSpqpwqKiqUZZ988gnS09ORkZGBEydOwGQy4bnnnsONGzecGPH9O3HihCrn/Px8AMCSJUuUdbRQ65s3byI0NBQZGRndLu9PfdevX4/s7GxkZWWhpKQETU1NmD9/Ptra2hyVxoD1lvetW7dQVlaGlJQUlJWV4cCBAzh37hwWLlzYZd1Vq1ap5sCOHTscEf4966veQN/zWmv1BqDKt6amBrt27YJOp8MLL7ygWs+V6t2fc5ZDj29xcU888YSsXr1aNRYcHCybN292UkSDr66uTgBIcXGxMpaYmChxcXHOC2oQpKamSmhoaLfL2tvbxWQyybZt25Sx5uZmMRqN8tVXXzkoQsdYt26dBAUFSXt7u4hos9YAJDs7W/m9P/W9fv266PV6ycrKUta5cuWKjBgxQg4dOuSw2O9H57y78+uvvwoAuXjxojI2e/ZsWbdu3eAGN4i6y7uveT1c6h0XFydRUVGqMVevd+dzlqOPb5e+otPa2orS0lJER0erxqOjo3H06FEnRTX4GhoaAKDLC9CKioowYcIETJkyBatWrUJdXZ0zwrOrqqoqmM1mBAYGIiEhARcuXAAAWK1W1NbWqmpvMBgwe/ZsTdW+tbUV+/btwyuvvAKdTqeMa7HWd+tPfUtLS2Gz2VTrmM1mWCwWTc2BhoYG6HQ6PPTQQ6rxb775BuPGjcPUqVPx1ltvufyVTKD3eT0c6n316lXk5eXh1Vdf7bLMlevd+Zzl6OPbpV/q+ffff6OtrQ0+Pj6qcR8fH9TW1jopqsElIti4cSMiIyNhsViU8djYWCxZsgQBAQGwWq1ISUlBVFQUSktLXfYpmzNnzsTevXsxZcoUXL16FR9++CEiIiJQWVmp1Le72l+8eNEZ4Q6KnJwcXL9+HStXrlTGtFjrzvpT39raWri7u2PMmDFd1tHK8d/c3IzNmzfjpZdeUr3kcdmyZQgMDITJZMLp06exZcsWnDp1SvlvTlfU17weDvXes2cPvLy8sHjxYtW4K9e7u3OWo49vl250Otz9L13gvz9s5zGtSEpKwu+//46SkhLVeHx8vPKzxWLBjBkzEBAQgLy8vC4HjauIjY1Vfp42bRrCw8MRFBSEPXv2KF9S1Hrtd+7cidjYWJjNZmVMi7Xuyb3UVytzwGazISEhAe3t7cjMzFQtW7VqlfKzxWLB5MmTMWPGDJSVlWH69OmODtUu7nVea6XeALBr1y4sW7YMI0eOVI27cr17OmcBjju+Xfq/rsaNGwc3N7cu3V1dXV2XTlEL1q5di9zcXBQWFsLPz6/XdX19fREQEICqqioHRTf4PD09MW3aNFRVVSl3X2m59hcvXkRBQQFee+21XtfTYq37U1+TyYTW1lb8+++/Pa7jqmw2G5YuXQqr1Yr8/HzV1ZzuTJ8+HXq9XlNzoPO81nK9AeDIkSM4e/Zsn8c74Dr17umc5ejj26UbHXd3d4SFhXW5fJefn4+IiAgnRWV/IoKkpCQcOHAAP//8MwIDA/vcpr6+HpcuXYKvr68DInSMlpYWnDlzBr6+vspl3Ltr39raiuLiYs3Ufvfu3ZgwYQLmzZvX63parHV/6hsWFga9Xq9ap6amBqdPn3bpOdDR5FRVVaGgoABjx47tc5vKykrYbDZNzYHO81qr9e6wc+dOhIWFITQ0tM91h3q9+zpnOfz4vtdvUQ8VWVlZotfrZefOnfLHH3/I+vXrxdPTU/78809nh2Y3r7/+uhiNRikqKpKamhrlc+vWLRERuXHjhmzatEmOHj0qVqtVCgsLJTw8XCZOnCiNjY1Ojv7ebdq0SYqKiuTChQty/PhxmT9/vnh5eSm13bZtmxiNRjlw4IBUVFTIiy++KL6+vi6dc4e2tjbx9/eX5ORk1biWan3jxg0pLy+X8vJyASDp6elSXl6u3F3Un/quXr1a/Pz8pKCgQMrKyiQqKkpCQ0Plzp07zkqrT73lbbPZZOHCheLn5ycnT55UHe8tLS0iInL+/Hl5//335cSJE2K1WiUvL0+Cg4Plsccec9m8+zuvtVbvDg0NDTJq1CjZvn17l+1dsd59nbNEHHt8u3yjIyLy5ZdfSkBAgLi7u8v06dNVt11rAYBuP7t37xYRkVu3bkl0dLSMHz9e9Hq9+Pv7S2JiolRXVzs38PsUHx8vvr6+otfrxWw2y+LFi6WyslJZ3t7eLqmpqWIymcRgMMisWbOkoqLCiRHbz48//igA5OzZs6pxLdW6sLCw23mdmJgoIv2r7+3btyUpKUm8vb3Fw8ND5s+fP+T/Fr3lbbVaezzeCwsLRUSkurpaZs2aJd7e3uLu7i5BQUHy5ptvSn19vXMT60Nvefd3Xmut3h127NghHh4ecv369S7bu2K9+zpniTj2+Nb9X1BEREREmuPS39EhIiIi6g0bHSIiItIsNjpERESkWWx0iIiISLPY6BAREZFmsdEhIiIizWKjQ0RERJrFRoeIiIg0i40OERERaRYbHSIiItIsNjpERESkWWx0iIiISLP+B1ta/5+zuv2BAAAAAElFTkSuQmCC", "text/plain": [ "
" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "for i in range(5):\n", " integrator.step(20000)\n", " plot.plot(integrator.getAdaptedFriction(0), label=f'{i}')\n", "plot.xlim([0, 200])\n", "plot.legend()\n", "plot.show()" ] }, { "cell_type": "markdown", "id": "9168045d-4048-49b3-b38a-d5e59bdb0f95", "metadata": {}, "source": [ "We can see the spectrum is well converged, so now we compute the RDF as before." ] }, { "cell_type": "code", "execution_count": 13, "id": "5ff1d432-b0b1-4869-9ef9-fcadd17e8863", "metadata": {}, "outputs": [], "source": [ "qtb_rdf = compute_rdf(context)" ] }, { "cell_type": "markdown", "id": "6ccfc38e-f4f3-4c5e-b3ff-489fed158dec", "metadata": {}, "source": [ "Now that we have all three versions of the RDF, let's compare them." ] }, { "cell_type": "code", "execution_count": 14, "id": "c71881a6-1386-49d1-bff3-4dc1ca197d0a", "metadata": {}, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAiMAAAGdCAYAAADAAnMpAAAAOnRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjEwLjMsIGh0dHBzOi8vbWF0cGxvdGxpYi5vcmcvZiW1igAAAAlwSFlzAAAPYQAAD2EBqD+naQAAZ8hJREFUeJzt3Xd4VHXaxvHvmZn0TkIaBAi9hCYoTZooiMra1rK66tpW7IoVe9tFd11fZBV7WdeuoOsiIqxUpZfQO4FQEtIgvc3Mef+YTCDUJEwyKffnuuYyOefMzDMHSG5/1TBN00RERETESyzeLkBERESaN4URERER8SqFEREREfEqhRERERHxKoURERER8SqFEREREfEqhRERERHxKoURERER8SqbtwuoDqfTyYEDBwgJCcEwDG+XIyIiItVgmib5+fnEx8djsZy8/aNRhJEDBw6QkJDg7TJERESkFvbu3Uvr1q1Per5RhJGQkBDA9WFCQ0O9XI2IiIhUR15eHgkJCZW/x0+mUYQRd9dMaGiowoiIiEgjc7ohFhrAKiIiIl6lMCIiIiJepTAiIiIiXtUoxoyIiEjzYJomdrsdh8Ph7VKkGqxWKzab7YyX3VAYERGRBqGsrIy0tDSKioq8XYrUQGBgIHFxcfj6+tb6NRRGRETE65xOJykpKVitVuLj4/H19dUilw2caZqUlZWRmZlJSkoKnTp1OuXCZqeiMCIiIl5XVlaG0+kkISGBwMBAb5cj1RQQEICPjw979uyhrKwMf3//Wr2OBrCKiEiDUdv/sxbv8cSfmf7URURExKsURkRERMSrFEZERETqgWEYfP/993X+PvPnz8cwDA4fPuyR19u9ezeGYZCcnOyR1zsRhREREREPSE9P595776V9+/b4+fmRkJDAuHHj+OWXX+q1jsGDB5OWlkZYWFi9vu+Z0GwaEWnWftl8kOyCMq7q31pTSaXWdu/ezZAhQwgPD+dvf/sbvXr1ory8nJ9//pm7776bLVu21Fstvr6+xMbG1tv7eYJaRkSk2TpUWMb4T1fx6LR1TJy+HofT9HZJchTTNCkqs3vlYZo1+7tw1113YRgGy5cv5/e//z2dO3emR48eTJgwgaVLl57wOY899hidO3cmMDCQ9u3b8/TTT1NeXl55fu3atYwcOZKQkBBCQ0Pp168fK1euBGDPnj2MGzeOiIgIgoKC6NGjBzNnzgRO3E3z22+/MXz4cAIDA4mIiGDMmDEcOnQIgFmzZnHuuecSHh5OZGQkl1xyCTt37qzR5z9TahkRkWZr9qZ0yh2uXzpfrtjL4aJyXv9DH/xsVi9XJgDF5Q66P/OzV9570wtjCPSt3q/InJwcZs2axV/+8heCgoKOOx8eHn7C54WEhPDxxx8THx/P+vXruf322wkJCeHRRx8F4Prrr6dv37689dZbWK1WkpOT8fHxAeDuu++mrKyMhQsXEhQUxKZNmwgODj7h+yQnJzNq1ChuueUWpkyZgs1mY968eZVL7hcWFjJhwgR69uxJYWEhzzzzDJdffjnJycn1NtVaYUREmq0Z69IAGNa5JUt3ZjNrYzq3fLyCd27oT7CffjxK9ezYsQPTNOnatWuNnvfUU09Vft2uXTseeughvvrqq8owkpqayiOPPFL5up06daq8PjU1lSuvvJKePXsC0L59+5O+z9/+9jf69+/P1KlTK4/16NGj8usrr7yyyvUffPAB0dHRbNq0iaSkpBp9ptrSvzYRaZZyCstYvDMbgOd/14MDh4v58ycr+W1HNte9t5SPbz6HFkG132tDzlyAj5VNL4zx2ntXl7tLp6Zjjr799lsmT57Mjh07KCgowG63ExoaWnl+woQJ3Hbbbfz73//m/PPP56qrrqJDhw4A3Hfffdx5553Mnj2b888/nyuvvJJevXqd8H2Sk5O56qqrTlrHzp07efrpp1m6dClZWVk4nU7AFXjqK4xozIiINEs/b0zH4TTpER9KYlQQQzpG8fntA4kI9GHdvlyuensxBw4Xe7vMZs0wDAJ9bV551CRYdOrUCcMw2Lx5c7Wfs3TpUq699lrGjh3LjBkzWLNmDU8++SRlZWWV1zz33HNs3LiRiy++mLlz59K9e3e+++47AG677TZ27drFDTfcwPr16+nfvz///Oc/T/heAQEBp6xl3LhxZGdn895777Fs2TKWLVsGUKWWuqYwIiLN0sz1ri6ai3rGVR7rnRDON+MHERfmz87MQn7/1mJ2ZhZ4q0RpJFq0aMGYMWN48803KSwsPO78idb7+O2332jbti1PPvkk/fv3p1OnTuzZs+e46zp37syDDz7I7NmzueKKK/joo48qzyUkJDB+/HimT5/OQw89xHvvvXfC+nr16nXS6cXZ2dls3ryZp556ilGjRtGtW7fKga31SWFERJqd7ILSyi6ai48KIwAdo0P49s7BtG8ZxIHcEq56ewnr9+V6o0xpRKZOnYrD4eCcc85h2rRpbN++nc2bNzNlyhQGDRp03PUdO3YkNTWVL7/8kp07dzJlypTKVg+A4uJi7rnnHubPn8+ePXv47bffWLFiBd26dQPggQce4OeffyYlJYXVq1czd+7cynPHmjhxIitWrOCuu+5i3bp1bNmyhbfeeousrCwiIiKIjIzk3XffZceOHcydO5cJEybUzU06BYUREWl2ft54EIfTJKlVKO2ijp/90Co8gG/uGETPVmHkFJZx7btLWLwzywuVSmORmJjI6tWrGTlyJA899BBJSUlccMEF/PLLL7z11lvHXX/ppZfy4IMPcs8999CnTx8WL17M008/XXnearWSnZ3NjTfeSOfOnbn66qsZO3Yszz//PAAOh4O7776bbt26ceGFF9KlS5cqA1SP1rlzZ2bPns3atWs555xzGDRoEP/5z3+w2WxYLBa+/PJLVq1aRVJSEg8++CB///vf6+YmnYJh1nQytRfk5eURFhZGbm5ulcE9IiK18cf3l/HrjiwevbALd43oeNLrCkrt3P6vlSzZlY2vzcJP9w+lQ8sTT5+UM1NSUkJKSgqJiYm13oZevONUf3bV/f2tlhERaVZcXTSuVo5ju2iOFexn46Obz6Z/2wjK7E5+rJgKLCKepTAiIs3KrI3pOE3o2SqMtpHHd9Ecy9/HyhVntQZg4bbMui5PpFlSGBGRZsU9i+biXqduFTnasM5RAKxOPURuUflprhaRmlIYEZFmI6uglCUnmUVzKq0jAunQMginCb9pIKuIxymMiEizMWuDq4umV+swEloE1ui5wztHA+qqEakLCiMi0mxUdtHUoFXEzd1Vs2BbZo13dBWRU1MYEZFmITO/lKW7XF00F9UijAxsH4mfzUJabgnbM7Qqq4gnKYyISLPgnkXTuxZdNOCaVTOgfSSgrhoRT1MYEZFmYea6ms+iOdawTke6akTEcxRGRKTJy8gvYVmKq4tmbFLtw8iILi0BWLYrh6Iyu0dqk8bvT3/6E4ZhYBgGNpuNNm3acOedd1bZcK5du3aV1wQGBpKUlMQ777xTef7jjz/GMIwT7i/z9ddfYxgG7dq1O+56wzCwWq1EREQwYMAAXnjhBXJzG99eSgojItLk/Vwxi6Z3QnitumjcOrQMplV4AGUOJ8t25XiwQmnsLrzwQtLS0ti9ezfvv/8+//3vf7nrrruqXPPCCy+QlpbGunXruOyyyxg/fjxfffVV5fmgoCAyMjJYsmRJled9+OGHtGnT5rj3DA0NJS0tjX379rF48WL+/Oc/88knn9CnTx8OHDhQNx+0jiiMiEiT92PFLJpLajFw9WiGYTCss6t1RF01cjQ/Pz9iY2Np3bo1o0eP5pprrmH27NlVrgkJCSE2NpaOHTvy0ksv0alTJ77//vvK8zabjeuuu44PP/yw8ti+ffuYP38+11133XHvaRgGsbGxxMXF0a1bN2699VYWL15MQUEBjz76aJ191rqgMCIiTZqri8bVijG2Z+wZv97wiim+GsRaD0wTygq98ziD6du7du1i1qxZ+Pj4nPI6f39/ysurruh766238tVXX1FUVAS4umMuvPBCYmJiqvXe0dHRXH/99fzwww84HI7afQAvsNXk4kmTJjF9+nS2bNlCQEAAgwcP5pVXXqFLly6nfN6CBQuYMGECGzduJD4+nkcffZTx48efUeEiItUxa0M6pgl9EsJpHVH7Lhq3wR2jsFoMdmUVsjen6Iy6feQ0yovgr/Heee8nDoDv6fcucpsxYwbBwcE4HA5KSkoAeO211054rd1u59NPP2X9+vXceeedVc716dOHDh068O2333LDDTfw8ccf89prr7Fr165q19K1a1fy8/PJzs4mOjq62s/zphq1jCxYsIC7776bpUuXMmfOHOx2O6NHj6awsPCkz0lJSeGiiy5i6NChrFmzhieeeIL77ruPadOmnXHxIiKnM6NiFs0lZzCL5mih/j70axMBqKtGjhg5ciTJycksW7aMe++9lzFjxnDvvfdWueaxxx4jODiYgIAA7r77bh555BHuuOOO417rlltu4aOPPmLBggUUFBRw0UUX1agW96J8hmHU/gPVsxq1jMyaNavK9x999BHR0dGsWrWKYcOGnfA5b7/9Nm3atGHy5MkAdOvWjZUrV/Lqq69y5ZVX1q5qEZFqyMgrYcVudxeNZ8IIuFZjXb47hwXbMvnjwLYee105hk+gq4XCW+9dA0FBQXTs2BGAKVOmMHLkSJ5//nlefPHFymseeeQR/vSnPxEYGEhcXNxJw8L111/Po48+ynPPPceNN96IzVajX9Vs3ryZ0NBQIiMja/Q8b6rZJzyGe/pQixYtTnrNkiVLGD16dJVjY8aM4YMPPqC8vPyEfWqlpaWUlpZWfp+Xl3cmZYpIM/VTRRdN3zbhtAoP8NjrDu8czauzt7F4RxZldie+Ng2/qxOGUaOukobk2WefZezYsdx5553Ex7u6mqKioioDy6m0aNGC3/3ud3z99de8/fbbNXrfjIwMPv/8cy677DIslsbz97LWlZqmyYQJEzj33HNJSko66XXp6enHDbyJiYnBbreTlXXi3S8nTZpEWFhY5SMhIaG2ZYpIM/bjutrvRXMqPeJDiQzypbDMwerUQ6d/gjQ7I0aMoEePHvz1r3+t1fM//vhjsrKy6Nq160mvMU2T9PR00tLS2Lx5Mx9++CGDBw8mLCyMl19+ubale0Wtw8g999zDunXr+OKLL0577bFNUafrz5o4cSK5ubmVj71799a2TBFppjLzS1mxx9VFU5u9aE7FYtEUXzm9CRMm8N5779Xqd1hAQMBpu1ny8vKIi4ujVatWDBo0iHfeeYebbrqJNWvWEBfn2b/zdc0wa7H95L333sv333/PwoULSUxMPOW1w4YNo2/fvrz++uuVx7777juuvvpqioqKTjv1CVw3PCwsjNzcXEJDQ2tarog0Q4u2Z3LDB8vp0DKIXx4a4fHX/27NPh78ai3d40KZef9Qj79+c1NSUkJKSgqJiYn4+/t7uxypgVP92VX393eNWkZM0+See+5h+vTpzJ0797RBBGDQoEHMmTOnyrHZs2fTv3//agUREZHaSM1xrdPQLrJuxhwM7eRqGdmUlkdGfkmdvIdIc1GjMHL33Xfz6aef8vnnnxMSEkJ6ejrp6ekUFxdXXjNx4kRuvPHGyu/Hjx/Pnj17mDBhQmWf1gcffMDDDz/suU8hInKMvTmun0t1tQ5IVLAfPVuFAbBo24nHv4lI9dQojLz11lvk5uYyYsQI4uLiKh9Hr62flpZGampq5feJiYnMnDmT+fPn06dPH1588UWmTJmiab0iUqf2VrSM1OWiZMM6axdfEU+o0dTe6gwv+fjjj487Nnz4cFavXl2TtxIROSPubpo2dRhGhneO5s15O1m0PROH08RqaTyLTIk0JI1nErKISA3URxjp2yacED8bh4rK2bC/8W3bLtJQKIyISJOTW1xObrFrA7LWEZ5b7OxYPlYLgzu6pl+qq0ak9hRGRKTJcY8XiQr2JcjvjBaaPq3hnV0bkSmMiNSewoiINDn1MXjVzT2IdU3qIXKLyk9ztYiciMKIiDQ59TFexK11RCAdo4NxmvDbTk3xFakNhRERaXL2HqpoGYmo+zACMKxiAbQFW9VVI1IbCiMi0uSkVix4Vh8tIwDDu7jCyMLtmdVaAkFk9+7dGIZBcnKyt0tpEBRGRKTJqc8xIwADElvgZ7OQllvC9oyCenlPaZr+9a9/cc455xAUFERISAjDhg1jxowZlef/9Kc/YRjGKR8nui4yMpILL7yQdevWeeujnZLCiIg0KQ6nyf5DFS0jkfUTRvx9rAxo75riO39rRr28pzQ9Dz/8MHfccQdXX301a9euZfny5QwdOpRLL72UN954A4DXX3+dtLS0ygfARx99dNwxgAsvvLDy2C+//ILNZuOSSy7xymc7nbqd8yYiUs8O5pVQ5nDiYzWIDa2/3V9HdG7Jwm2ZzNuSyZ+Hdai3923KTNOk2F58+gvrQIAtoLKVoTpmzZrFSy+9xIYNG7BarQwaNIjXX3+dDh1cfxeWL1/OHXfcwebNm0lKSuLJJ5+s8vylS5fyj3/8gylTpnDvvfdWHv/LX/5CSUkJEyZM4NJLLyUhIYGwsLAqzw0PDyc2Nva4mvz8/CqPx8bG8thjjzFs2DAyMzNp2bJltT9bfVAYEZEmxT2TplV4QL0uzz6qWzQvzNjEit055JWUE+qvXcnPVLG9mAGfD/DKey+7bhmBPtVvWSssLGTChAn07NmTwsJCnnnmGS6//HKSk5MpLi7mkksu4bzzzuPTTz8lJSWF+++/v8rzv/jiC4KDg7njjjuOe+2HHnqI1157jWnTpvHAAw/U6vMUFBTw2Wef0bFjRyIjI2v1GnVJYUREmpTUeh4v4tY2Moj2LYPYlVnIom1ZXNwrrl7fX7zr2M1fP/jgA6Kjo9m0aROLFy/G4XDw4YcfEhgYSI8ePdi3bx933nln5fXbtm2jQ4cO+Pr6Hvfa8fHxhIWFsW3bthrVNGPGDIKDgwFXWIqLi2PGjBlYLA1vhIbCiIg0Kfu8FEYARnWNZldmCr9sOagw4gEBtgCWXbfMa+9dEzt37uTpp59m6dKlZGVl4XQ6AUhNTWXz5s307t2bwMAjfycHDRpUo9c3TfOEQeVURo4cyVtvvQVATk4OU6dOZezYsSxfvpy2bdvW6LXqmsKIiDQp9bng2bFGdo3mvUUpLNiqXXw9wTCMGnWVeNO4ceNISEjgvffeIz4+HqfTSVJSEmVlZdWa7t2pUyd+/fVXysrKjgsdBw4cIC8vj86dO9eopqCgIDp27Fj5fb9+/QgLC+O9997jpZdeqtFr1bWG11YjInIGvBlGzm7XghA/G9mFZazdd7je31+8Izs7m82bN/PUU08xatQounXrxqFDhyrPd+/enbVr11JcfGQw7tKlS6u8xh/+8AcKCgp45513jnv9V199FX9/f6655pozqtMwDCwWS5U6GgqFERFpUup7wbOj+VgtDOvsmqUwb4um+DYXERERREZG8u6777Jjxw7mzp3LhAkTKs9fd911WCwWbr31VjZt2sTMmTN59dVXq7zGoEGDuP/++3nkkUf4xz/+wc6dO9myZQtPPfUUU6ZM4b333qvxwNPS0lLS09NJT09n8+bN3HvvvRQUFDBu3DiPfG5PUhgRkSajuMxBVkEpUH9LwR/rvK6uXXx/2aww0lxYLBa+/PJLVq1aRVJSEg8++CB///vfK88HBwfz3//+l02bNtG3b1+efPJJXnnlleNeZ/LkyUydOpUvvviCpKQkunXrxt///nfmzp3LH//4xxrXNWvWLOLi4oiLi2PAgAGsWLGCb775hhEjRpzJx60ThtkI1i7Oy8sjLCyM3NxcQkNDvV2OiDRQ2w7mM/r/FhLqb2Pdc2O8UkN2QSn9//I/TBOWThxFbFj9rXXSmJWUlJCSkkJiYiL+/rpn4Foyfvjw4QwaNIjPPvsMq9Xq7ZJO6FR/dtX9/a2WERFpMlKzK8aL1NPKqycSGexHn4RwAOaqq0bOQLt27Zg/fz5du3Zt8nvYKIyISJNRucaIl7po3M7r4uqqURiRM5WYmMhzzz1Hv379vF1KnVIYEZEmY+8h782kOdrIinEjv+3IoqTc4dVaRBoDhRERaTLqe7fek+kRH0pMqB/F5Q6W7sr2ai0ijYHCiIg0Gd5cY+RohmFUzqrRFN+aaQRzKuQYnvgzUxgRkSbBNE32Vqwx4u2WEYDzusYA8MuWDP2CrQYfH9fGgkVFRV6uRGrK/Wfm/jOsDS0HLyJNQlZBGcXlDgzDtWOvtw3pGImvzcK+Q8XsyCigU0yIt0tq0KxWK+Hh4WRkuFqSAgMDMQwtp9+QmaZJUVERGRkZhIeHn9HUY4UREWkS3F008WEB+Nq83+gb6GtjUPtIFmzL5JctGQoj1RAbGwtQGUikcQgPD6/8s6sthRERaRKODF71fquI23ldo1mwLZO5WzIYP7yDt8tp8AzDIC4ujujoaMrLy71djlSDj4+PRxZjUxgRkSZhbwNZY+Ro53WN5tkfNrJqzyFyi8oJC6x9n3pzYrVaG+xqo1I3vN+WKSLiAQ1lJs3REloE0ik6GIfTZMH2TG+XI9JgKYyISJNQGUa8uBT8iZzXrWI11s0HvVyJSMOlMCIiTYK7m6Z1A+qmgSNLw8/flonDqSm+IieiMCIijV6Z3UlaXgnQsLppAPq1jSDU38bhonLWpB7ydjkiDZLCiIg0evsPF2OaEOBjJSrY19vlVGGzWhiujfNETklhREQavaMHrzbEhbLO69oSUBgRORmFERFp9FIb4BojRxveORqLAVvS89l/uNjb5Yg0OAojItLo7Wsgu/WeTIsgX/q2iQDUOiJyIgojItLoNcQ1Ro6lXXxFTk5hREQavdQGuPrqsdxh5LcdWRSXObxcjUjDojAiIo1eQ13w7GhdY0OID/On1O5kya4sb5cj0qAojIhIo5ZbVE5+iR1o2C0jhmEwsqJ15JfN6qoROZrCiIg0au5WkZYhfgT4NuzN1dxdNfO3ZmKaWo1VxE1hREQatSPjRRrmtN6jDeoQia/Vwv7DxezIKPB2OSINhsKIiDRqjWEmjVugr40B7VsArtYREXFRGBGRRm3vocYTRgCGd3atxjp/m8aNiLgpjIhIo7a3gS94dqwRFfvUrEg5RGGp3cvViDQMCiMi0qilNrIw0qFlEK0jAihzOFm8M9vb5Yg0CAojItJoOZwm+w+59nppLN00hmEwoktFV81WddWIgMKIiDRiabnF2J0mvlYLMaH+3i6n2kZ20RRfkaMpjIhIo+XuomkVEYDVYni5muo7eorvzkxN8RVRGBGRRquxDV510xRfkaoURkSk0dqb4x4v0vAXPDtW5RRfhRERhRERabwa04Jnx3JP8V2ekqMpvtLsKYyISKN1ZCn4xhdGjp7iu0RTfKWZUxgRkUarQY0ZKcqBjd/DgWRwOk57eZUpvlqNVZo5m7cLEBGpjcJSO9mFZQC0ifRyGNnxC+b3d5FXdJAwpwm+IZBwDrQdBG0GQ6t+4HP81OMRnaP5dGlq5RRfw2g8M4JEPElhREQaJfeeNGEBPoT6+3iniPJiHHOeZe6GT/ggPJSNLRMYW1TGxMwMInb+Ajt/cV1n9YX4s46Ek8Rh4OPP4I6uKb77DhWzM7OQjtHB3vkcIl6mMCIijVJqtncHr5bvX82MGbfzoaWQ3TEtK4//FOjLsk7debrlEM7POQipS6DgIOxd6nrwfxDaGs57ksBe13BOYgt+3ZHF/K0ZCiPSbGnMiIg0St6aSVNUWsC/ZtzGhbP+yDP+Zez29SHEGsCfe/2Zdy54hw5hHcgpy+PB/T/xSEw0h+5eCveuhkvfhD5/hOBYyNsH398Jb5/LDZFbAFNTfKVZU8uIiDRK+yr2pKmvwauHSg7xefLbfL7lS/IMJ9isRGPjxl638/ukmwjyCQLg63Ff8/bat/lww4fM2j2L5enLeWrgU1zQ94/Q949QXgzL34VF/4CMTYzJuI+vfLvyj93XUVjajyA//ViW5kctIyLSKB3ZrbfuFzz7bPNnjP76PN7e+jl5hpO25Q6ejx/DT9cv5aa+d1UGEQBfqy/3nXUfn130GR3DO5JTksOE+RN4eMHD5JTkgE8ADLkf7l8LQx7AtPkzwLKFr23PUPjJtZC5rc4/j0hDozAiIo1SfXXT/G/P/3h5+cuUmHa6l5byj/JQ/nPZ91xxwav42vxO+rweUT346pKvuL3n7VgNKz/v/pnL/3M5s3fPdl0QEAEXPI9x72pWtrgEh2kQvX8OTB0AP9wL+el1+rlEGhKFERFpdEzTrFxjpC7DyNacrTyx6HEArs/N58sONzL65vlYozpX6/mVrSQXH2kleWjBQ0xcNJGiclf9hLXi8PmvMabsFRZazgHTCas/gffPh+JDdfXRRBoUhRERaXQy80sptTuxGBAfXjfdNDklOdw39z6KHaUMLC7m4bBeGKOeAmvNpxH3iHS1ktzR6w6shpUZu2Zw/czrSclNAWBwx0hSLW24segB9l3+PYS3hdy9MONBME0PfzKRhkdhREQaHXcXTVxYAD5Wz/8YK3eUM2H+BA4UHiCh3M6rGdnYRkw8o9f0tfpyT997eH/0+0QFRLHj8A6unXEts3bPItDXxjmJrl18Z+W1hd9/BIYVNn4Ha7/0xEcSadAURkSk0XEveFZXg1dfXv4yqw6uIsiw8s+DGYS1GeJasMwD+sf255tx33B27NkU2Yt4ZMEjvLL8FYZ2DgdgwbZMaN0PRlaEn5kPQ06KR95bpKFSGBGRRicjrxRwtYx42ldbvuLrbV9jYPDKwSw6lNth+GMefY+ogCjeveBdbk26FYBPN3/KnEPPYdhyWbYrh6IyO5w7AdoMgrICmP5ncGhnX2m6FEZEpNFx70kTGeTr0dddnracl5e/DMD9QZ0YXpgPbYdA4lCPvg+AzWLjgX4PMGXkFEJ8Qth6eAMh7afg8Nvq2sXXYoUr3gW/UNi3HBa96vEaRBqKGoeRhQsXMm7cOOLj4zEMg++///6U18+fPx/DMI57bNmypbY1i0gzl1XgahmJDD751Nqa2pu/l4cWPITdtHNRwnncsnmR68TwRz32Hicyss1Ivhr3Fd1adMO0FhLQ5kPeXvsOTtMJ4W3g4tdcFy54BVKX1WktIt5S4zBSWFhI7969eeONN2r0vK1bt5KWllb56NSpU03fWkQEgOyCipaRYM+0jBSWF3Lf3Ps4XHqY7pHdeb7UH8NRCgkDIXG4R97jVBJCEvj3Rf9mYMuLMAyTLaXf8OSvT7pO9roKel7tmvI7/XYoyavzekTqW43DyNixY3nppZe44ooravS86OhoYmNjKx9Wq7Wmby0iAkB2oatlJMoDYcRpOnli0RPsOLyDqIAoXj/nGfxX/ct1cvijYBhn/B7V4Wf1Y/Kov1CefjWmaWHGrhlsyt7kOnnxqxDWBg7vgZ88O35FpCGotzEjffv2JS4ujlGjRjFv3rxTXltaWkpeXl6Vh4iIW2XLSNCZd9O8tfYt5u6di4/Fh8kjJxOb/AXYi6FVf+hw3hm/fk0E+dnoF3kB9rxeAHy84WPXCf8w1/gRwwJrP4cN0+q1LpG6VudhJC4ujnfffZdp06Yxffp0unTpwqhRo1i4cOFJnzNp0iTCwsIqHwkJCXVdpog0EqZpHhnAeoYtI+mF6by37j0Anh30LL0DW8GKD1wnRzxeb60iRxvROZqy7GEA/LznZ/bl73OdaDsIhj7s+nrGg3B4b73XJlJX6jyMdOnShdtvv52zzjqLQYMGMXXqVC6++GJeffXkI8MnTpxIbm5u5WPvXv2jExGXglI7ZXYncOYtI19t/QqH6eDs2LO5tOOlsPifUF4E8X2h4/meKLfGRnaNxlkaj6OwM07Tyb82/uvIyeGPQqt+UJIL340Hp8MrNYp4mlem9g4cOJDt27ef9Lyfnx+hoaFVHiIicKSLJtDXSoBv7ceelTpKmbbN1d1xXdfroDAblrtaSRjunVYRgI7RwXSOCaY0yzVw9vsd37t2+wXXUvRXvAc+QbDnV1g8xSs1iniaV8LImjVriIuL88Zbi0gj5x68eqZdNLNSZnGo9BCxQbGMSBgBS9+E8kKI7QWdx3ig0tq7uGc8jqL2BJrtKHGU8MWWL46cjOwAF/3N9fXclyB7p3eKFPGgGoeRgoICkpOTSU5OBiAlJYXk5GRSU1MBVxfLjTfeWHn95MmT+f7779m+fTsbN25k4sSJTJs2jXvuucczn0BEmpUsDwxeNU2Tz7d8DsA1Xa7BVpIHy951nRz+mNdaRdwu7hULGBw6MASAL7Z8cWSXX4A+10OHUeC0w/yXvVOkiAfVOIysXLmSvn370rdvXwAmTJhA3759eeaZZwBIS0urDCYAZWVlPPzww/Tq1YuhQ4fy66+/8uOPP9Z4arCICBzppjmTab1rM9eyKXsTvhZfrux0JSx7G8ryIaYndL3YU6XWWsfoELrGhlCW14MI3zhyS3P5bsd3Ry4wDBjl+pnL+m/g4CbvFCriIbaaPmHEiBGYp9jS+uOPP67y/aOPPsqjj9btCoYi0nxku1dfPYOWEXeryNjEsUSYBix923Vi+CNebxVxu7hnHFvS8wkuOZ9Dln/zycZPuLrL1fhYfFwXxPeBbr+DzT/AvL/AtZ95tV6RM6G9aUSkUTnTab2ZRZnM2T0HgOu6XQfL3oHSXIjuDl3HeazOM3VRL9e4um07uhDh14IDhQeYvXt21YtGPgkYsGUG7F9d/0WKeIjCiIg0Kme6L803277Bbtrp07IP3UPawdKprhPDHgFLw/mR2KFlMN3iQrE7bPQOdXUdfbTho6ot09Fdodc1rq/nvuSFKkU8o+H8yxMRqYYzGTNS7ijnm23fABWtIhunQ8lhiGgH3S/1YJWecUlF60h2en8CbAFsPbSVxQcWV71oxONgscHOX2DP4hO8ikjDpzAiIo1K5dTeWowZmb1nNlnFWbQMaMn5bc+HlR+6TvS7GSwNb7+si3q6wsiyHSWMS3QN+v9ow0dVL2qRCH1vcH39y4twijF9Ig2VwoiINCpnsmOve72Oq7pchU/6Rti/Ciw+0PePHq3RUxKjgugRH4rDaRJtno/NsLEsfRkbszZWvXDYI2D1g9TFrhYSkUZGYUREGg2H0+RQUe3CyMbsjazNXIvNYuOqzlcdaRXpfikERXm6VI+5uKKr5tctDsYmjgXgww0fVr0orBWcfZvr67kvqXVEGh2FERFpNA4XleGs+D0bEVizMPL5Ztd03tFtRxNl+MD6b10nzr7VkyV63MUVXTWLd2ZzefvrAfhf6v9IzUuteuG5D7qWiT+wxjW7RqQRURgRkUbDPa03PNAHH2v1f3zllOQwK2UWUDFwdd3XrqXfW3aFNoPqpFZPaRsZRFIrV1fN9n0hDG019PgN9ACCW8LA8a6v5/5Fm+hJo6IwIiKNRuW03qCatYpM2zaNMmcZPSJ70Cuy55Eumv63NJhFzk7l4p7xAPy4/gA3J90MuDbQyyrOqnrh4HvBLwwyN8OGafVdpkitKYyISKNxZPBq9WfS2J12vtr6FeBqFTH2LYeMTeATCL2vrZM6Pc3dVbNkZzbtgnrSK6oXZc6yqhvoAQREwJB7XV/P+ys4yuu5UpHaURgRkUbDvRR8TdYYmZs6l4NFB2nh34IL210IKz5wnUi6EvzD6qJMj2sTGUiv1mE4Tfh548HK1pEvtnzB3vy9VS8ecCcERsGhFEjWEvHSOCiMiEijUbkUfA3WGHG3HlzZ6Up8S/Jh0/euE/1v8XR5dcrdOvLjujRGJoykW4tu5Jflc/vs20kvTD9yoV8wDJ3g+nrB36C8xAvVitSMwoiINBpZNVxjZGvOVlYeXInVsHJ1l6sh+VNwlEF8X2h1Vl2W6nGVC6ClZJNdWM6bo96kTUgb9hfs57bZt5FZlHnk4v63Qkg85O2HVR+d5BVFGg6FERFpNLJruC+Nu1XkvDbnERsQDSsrfjE3slYRgIQWgfROCHd11WxIp2VgS94f/T7xQfHsydvD7bNvJ6ckx3Wxj79rB2KARf+A0gLvFS5SDQojItJouLtpoqoxmyavLI8fd/0IwHVdr4OU+a5xFH5hrvEijdAlFa0jM9alARAXHMf7Y94nOiCanbk7uWPOHeSW5rou7nuDa8+dwkxY8LKXKhapHoUREWk0atIysnDfQkocJbQPa0+/mH5HBq72vhZ8g+qyzDoztmcsAMt355CR5xoLkhCSwHtj3qOFfwu25Gzhzv/dSUFZAVh94MJXXE9c8qZrMTSRBkphREQajZrsS7Nw30IARiaMxMhPg60/uU70v7nO6qtrrSMC6dsmHNOEnzYcGbTaPqw9741+jzC/MNZnrefuX+6mqLwIulzoagUynfCfezXVVxoshRERaRRK7Q7yS+0ARJ1mNo3daefX/b8CMKz1MFj9CZgOaDMYorvVea116ehZNUfrHNGZdy54h2CfYFZnrOb+efdT6ih1tY4ERMDB9bD4n94oWeS0FEZEpFHIqRgvYrMYhAbYTnnt2sy15JflE+YXRq8W3WFVxdLpDXwfmupwz6pZsSeH9Nyq03Z7RPbgrfPfIsAWwNK0pUyYP4HygHAY81fXBfNfhuyd9VyxyOkpjIhIo+DuomkR5ItxmiXc3V00Q+KHYNv5C+QfgMBI6Dauzuusa/HhAfRrG1HRVZN23Pk+0X14c9Sb+Fn9WLhvIY8tegx7z6ug/UhwlMIP94HT6YXKRU5OYUREGoWsGg5ehYouGvfA1b5/BFv1F0tryE7WVeN2duzZvD7ydXwsPszZM4dPN38G4ya7lsDf8yus+aQeqxU5PYUREWkU3C0jp1sK/kDBAXYc3oHFsHBuYALs/MV1ol/jHbh6LHdXzco9h9h/uPiE1wxpNYQnBzwJwNS1Uzng4wvnPeU6OfsZyDtxkBHxBoUREWkUsgurt2Ovu1WkT8s+hK2v2Lm2wyhokVin9dWn2DB/BrZvAcD3a/af9LorOl3BWdFnUWwv5q/L/op5zh0QfxaU5sLMh+urXJHTUhgRkUahujv2usPI0PjBsOZT18FGuOLq6Vx5VmsApq3ah2maJ7zGMAyeGfQMNouNBfsW8Mu++fC7f4LFBltmwKYf6rFikZNTGBGRRqE6+9IU24tZnr4cgGFGEBRlQXAMdL6wXmqsT2N7xhHgY2VXViHJew+f9LoO4R24uYeri2rSskkUtGgHQx5wnZz5MBQfqvNaRU5HYUREGgV3N82p1hhZnracUkcpcUFxdDq4w3WwwyiwnnoqcGMU7GfjwiTXiqzTVu875bV/7vVnEkISyCjO4I3kN2DYIxDZEQoOwpxn6qNckVNSGBGRRqE6q68ePYvG2ONa9IzEoXVem7e4u2r+uzaNUrvjpNf52/x5aqBr8Ornmz9nY+5OGDfFdXL1J5CysM5rFTkVhRERaRROty+NaZos3F8RRmL6w/7VrhPtmm4YGdQhkthQf3KLy5m7OeOU1w6OH8xFiRdhYvL8kuextxlwZCzNf++H8hPPyhGpDwojItLgmaZZuWPvyWbTbDu0jfTCdPyt/pxT4nAt/x7RDsIT6rHS+mW1GFx+Vivg9F01AI+c/QghviFsztnMF1u+gPOfg5A4yNkF8yfVcbUiJ6cwIiINXmGZg1K7a9XQk3XTLNq/CIBz4s7Bf+8S18Em3CridmVFGJm/NbNyYbiTiQqI4sF+DwLwxpo3SHcUw8WvuU4ufkM7+4rXKIyISIPn7qIJ9LUS6HviwaiV40VaDYMUVzAhcVi91OdNHaND6N06DLvT5IfkA6e9/spOV9KnZR+K7EVMWjYJul4EPS53tSRpZ1/xEoUREWnwso7al+ZEDpccZm3mWgCGRfWB9HWuE82gZQTgyn4Va45Uo6vGYlhca48YNubuncvc1Lkw9m9H7ew7pa7LFTmOwoiINHinG7z664FfcZpOOkV0Ii5rJ5hOiOwEoXH1WabXjOsVj4/VYOOBPLak5532+k4Rnbipx00ATFo+iSK/YLjwZdfJ+a9A5ra6LFfkOAojItLguQevRp2kZaRqF03FNNUmPKX3WBFBvpzXNRqA6atPvjz80e7ofQetgluRXpjOm8lvQq9rXGuyOErhv9rZV+qXwoiINHhHWkaODyN2p53f9v8GVOzSu7tivEgz6aJxc6858t2a/dgdpw8SAbaAyo30Pt38KRuzN1Xs7BsEqUtg5Qd1Wa5IFQojItLgZZ1iX5q1mWvJK8sjzC+MXoGt4OAG14lmFkZGdImmRZAvmfmlLNqRVa3nDG09lAvbXYjTdPLA/AfI9A1wTfcF+N9zcHhvndUrcjSFERFp8E61xoi7i2ZI/BBse5e6DrbsBsEt662+hsDXZuF3veOB6nfVADw18CnahbYjvTCde+feS1Gf6yBhAJQVwIwH4SSb8Il4ksKIiDR47m6aqBO0jBy9BHxzmtJ7Iu6umtkb08krqd4U3TC/MKaOmkqEXwQbszfyxOKncI57Hay+sGMOrP+mLksWARRGRKQRONm+NAcKDrDj8A4shoVzW517ZLxIMxq8erSkVqF0jgmm1O7kx3Vp1X5eQmgCr5/3Oj4WH35J/YXJqT/B8EddJ396DAqr1+0jUlsKIyLS4Ll37I08Zsded6tIn5Z9CCsrgcwtgAFth9R3iQ2CYRhcUdE6Mr0aa44crW90X14a8hIAH238iG+j20BMEhTnuAKJSB1SGBGRBs3pNMlxT+09pmXEHUaGth56pFUkNgkCW9RrjQ3J5X1bYTFgxe5D7MkurNFzL2p/EXf1uQuAl5ZPYsmQ8WBYYMO3sHVWXZQrAiiMiEgDd7i4HGfFGMqIowawFtuLWZ6+HDh2Sm/zHC/iFhPqz7mdXIN3azKQ1W18r/Fc0v4SHKaDhza8xc7+N7pOzHgQSk6/oJpIbSiMiEiD5h68Gh7og4/1yI+sFekrKHWUEhcUR6fwTkcNXm2e40WO5t48b/qafTidNZsNYxgGzw9+nrOizyK/PJ+7izeT3aId5B+An5+og2pFFEZEpIE72b40C/YuAFytIkZ+GuTsdHUptB1c7zU2NKO7xxLsZ2NvTjErdufU+Pm+Vl9eH/k6bULasL8wjftat6XEsMCaf8M6za4Rz1MYEZEGzT14NeqowaumabJw/wmm9Mb1Af+w+i6xwQnwtXJxT9e+PLXpqgEI9w/nzVFvEuobyrr8FJ7uNhgnwIwHIGuHx2oVAYUREWngTjStd/vh7aQXpuNv9eec2HOa5X40p3NFRVfNj+vTKC5z1Oo12oW1Y/LIydgsNmYVp/Juu56uxdC++ROUF3uwWmnuFEZEpEE70b40Sw+4VlrtH9sff5s/7K4II8188OrRzm7XgoQWARSU2nl7wc7av07s2Twz8BkA3rEUsDU0Gg6uh1kTPVWqiMKIiDRsWZVLwR/pplmftR6As6LPgkN74HAqWGzQZqBXamyILBaDh0d3AeCfc7ezshZjR9wu73Q557c5H7vp4Lm2nXFgwKqPYP23nipXmjmFERFp0I4sBX+kZWRd5joAerXsdWRKb/xZ4Bdc7/U1ZJf2acUVfVvhNOH+L5PJLa7eEvEnMnHAREJ8QthQkMrnfce5Dv73fsiufauLiJvCiIg0aNnH7NibVZzFgcIDGBj0iOzR7PejOZ3nL+1BmxaB7D9czJPfrces5cZ30YHRTOg/AYB/5m9hX9uBrvEjX98E5SWeLFmaIYUREWnQjt2xd32mq4umQ3gHgn2Cmv1+NKcT4u/D69f2wWoxmLEujWm1nF0DcEWnK+gf059iRwkvxsZhBka5xo/8rPEjcmYURkSkQTsygNXVMrIu66gumpxdkLfftcNswgCv1djQ9W0TwYQLOgPwzH82sDurZsvEu1kMC88Nfg5fiy+LM1Yx49zbAQNWfggbpnmwYmluFEZEpMEqszvJK7EDR8aMuFtGekb1PDKlt/XZ4BPglRobi/HDOzAgsQVFZQ7u+3INZXZnrV6nbWhb7uxzJwCvpM4ge/DdrhM/aPyI1J7CiIg0WO4N8mwWg1B/HxxOBxuyNwAVYaRyPxp10ZyO1WLwf9f0ISzAh3X7cvm//22r9Wvd1OMmurboSm5pLq/4lbl2SS7Lh280fkRqR2FERBqsrIoumhZBvlgsBrtyd1FYXkiALYCOYR20H00NxYcH8PIVPQF4e8FOFu/IqtXr+Fh8eG7wc1gMCz/tnsXCIX+GwChI1/gRqR2FERFpsNyDV9370rjXF0mKSsKasxMKM8Dm7+qmkWoZ2zOOP5yTgGnCg18nc6jiHtdUj8ge3NjdtaPvi+vepPB3U6gcP5L8uQcrluZAYUREGqwja4xUDF6tWF+kyniRhAFg8zvh8+XEnr6kO+1bBnEwr5THpq2r9XTfu/rcRevg1qQXpvN67loY8bjrxIwH4UCy5wqWJk9hREQarGP3pakyk0ZTemst0NfGlGv74mM1mL3pIJ8vT63V6wTYAnh28LMAfLnlS5K7jobOF4K9BL66AYpqv+qrNC8KIyLSYGVV7NgbGeRHYXkhOw+7Zmv0apF0ZLyI9qOplaRWYTx2YVcAXpyxiR0Z+bV6nYFxA7ms42WYmDy79HnKfvdPaNEeclPh21vAWbtN+qR5URgRkQbr6JaRjVkbcZpOYoNiaVmQCcU54BMErc7ycpWN1y1DEhnaKYqScicTp9d+ddaH+z9MpH8ku3J38f6Ob+Gaz8AnEHbNg7kverhqaYoURkSkwTp6Xxp3F03PqJ6wZ7HrgjYDwOrjrfIaPYvF4JUrexHoa2XF7kO1Xp01zC+MiQNcs2jeXfcuy50FcOkbrpO//h9s+o+nSpYmSmFERBqs7KN27HUvdta7ZW/Yv8p1gVZdPWPx4QHcP6oTAJNmbuZwUe1m14xuO5rfdfgdDtPBIwsfIT1xCAy6x3Xy+7sgc6unSpYmSGFERBosdzdNiyCfqi0j7jDSqp+3SmtSbjk3kU7RwWQXlvH3n2sXGgzD4OmBT9OtRTdySnJ4cN6DlI58wrUgXVkBfHk9lOR5uHJpKhRGRKRBMk2T7IoBrKb1MFnFWVgNK92CWkH2dtdF8Rov4gk+VgsvXZYEwOfLU0nee7hWr+Nv8+e1Ea8R5hfGhuwNTFr5d/j9RxBa8Wf2/Z3grN0y9NK0KYyISINUVOagpNz1i+tAiev/1jtHdCYgY7Prgoh2EBTppeqangHtI7mibytME576fj0OZ+0Gs7YOac3fhv4NA4Np26fx7YEFcPW/XZsZbpkBv77m4cqlKVAYEZEGyd1FE+BjZcsh1340vVr2UhdNHZp4UTdC/G1s2J/HZ8v21Pp1BrcazH1n3QfAX5f9lfV+vnDRq66Tc1+CHf/zRLnShNQ4jCxcuJBx48YRHx+PYRh8//33p33OggUL6NevH/7+/rRv35633367NrWKSDNSucZIsG/lMvCu8SKrXRcojHhcyxA/Hh3TBYC//7yVjPzab3p3a9KtnJdwHuXOch6c/yDZ3S+Bs24CTPj2Vjhcu4XWpGmqcRgpLCykd+/evPHGG9W6PiUlhYsuuoihQ4eyZs0annjiCe677z6mTZtW42JFpPk4MnjVyqbsTYBaRurDdQPa0rNVGPkldibN3FLr1zEMg7+c+xfahbbjYNFBHl34KPYLJ7nG+ZQchm9uBnvtZu5I01PjMDJ27Fheeuklrrjiimpd//bbb9OmTRsmT55Mt27duO2227jlllt49dVXa1ysiDQf7jVGAoIzKHWUEuIbQlt8ID8NDCvE9vJyhU2T1WLw0mVJGAZ8t2Y/S3Zm1/q1gn2DmTxyMoG2QJanL+f1dW/DVR+BfxjsXwm/PO/ByqUxq/MxI0uWLGH06NFVjo0ZM4aVK1dSXl5+wueUlpaSl5dX5SEizYt7jRF8XWMXekb1xHJgjetYdHfwDfRSZU1f74Rwrh/QBoCn/7OBMnvtZ8B0CO/AS+e+BMDHGz9mVu4WuHSq6+SSN2DLzDOuVxq/Og8j6enpxMTEVDkWExOD3W4nKyvrhM+ZNGkSYWFhlY+EhIS6LlNEGpisipaRYutu4NguGk3prWuPjO5KZJAvOzIK+ODXlDN6rQvaXsDNSTcD8Mxvz7A9tgsMvMt18vs7NX5E6mc2jWEYVb53739w7HG3iRMnkpubW/nYu3dvndcoIg2Le8zIIYdrTREtdla/wgJ9eOKibgBM+WU7+w8Xn9Hr3df3PgbEDaDYXsy9c+8lufcVGj8ileo8jMTGxpKenl7lWEZGBjabjcjIE68R4OfnR2hoaJWHiDQv2YWlYCnicPkBAHq26AH7K7ppFEbqxRVnteKcdi0oLnfwwn83ntFr2Sw2/j7s77QKbsX+gv3cMPsWnu/cj9wAjR+ReggjgwYNYs6cOVWOzZ49m/79++Pjow2uROTEsgvKsAbsA6BNSBsiCrOgLN+1G2zLrl6urnkwDIMXL0vCZjH4eeNB5m45eEavF+EfwRcXf8FlHS8D4Ns9P/O71q34MSgQU+NHmrUah5GCggKSk5NJTk4GXFN3k5OTSU119flNnDiRG2+8sfL68ePHs2fPHiZMmMDmzZv58MMP+eCDD3j44Yc98wlEpEnKKijDGuD6udKz5VFdNHF9wGrzXmHNTJfYEG49NxGAidPXk1N4Zt0pEf4RvDjkRT4c8yGJYYnk2At4PDqKP8e2JPW/d2n8SDNV4zCycuVK+vbtS9++fQGYMGECffv25ZlnngEgLS2tMpgAJCYmMnPmTObPn0+fPn148cUXmTJlCldeeaWHPoKINDVOp8mhojKsAa7xYlXHi2jwan27//xOdGgZxMG8Uh75Zm3luL8zcXbs2Xw77lvu6XMPvhZflgYEcHlUMG9Pu5qy0gIPVC2NiWF64m9VHcvLyyMsLIzc3FyNHxFpBg4VltH3xdkEd3oRw1bEFxd/QdL0e+DAGtfGa0nVW+dIPGfTgTwum/obZXYnz47rzs1DEj322ql5qbz065MsyUwGoJ01mGfOn8LZsWd77D3EO6r7+1t704hIg5NdWIrhk41hK8LX4kuX4LaQ7tqfRoNXvaN7fChPVsyumTRzCxv253rstduEtuGdsZ/wSodriLQ72O0o4Jafb+HttW97pBVGGj6FERFpcFzjRVxdNF0ju+KTtRWc5RAYBeFtvFxd83XjoLZc0D2GMoeT+75YQ2Gp3WOvbRgGF537FD/EXcSV+a5umjeT32TirxMpdZR67H2kYVIYEZEGJ/uowau9oo7Zj+Yk6xNJ3TMMg79d2YvYUH92ZRXy3A9nNt33REJH/5Xn/NrzTFY2NtPkx10/cuvMG8gurv2y9HJ6+SUnXhG9viiMiEiDk11YWtkyos3xGpaIIF8mX9sHiwHfrNrHf5L3e/YNbL5ww3dc1f1G3jqYRYjDydqczVw3/RK2Z23y7HsJucXl/HXmZgZNmsvenCKv1aEwIiINTnp+ARb/NEArrzZEA9tHcs95nQB48rsN7Mku9Owb+IfC2JcZeNP/+IxY2pSXc8BewA0zrmHhiurtGC+nVu5w8q/Fuxnx93m8u3AXBaV2zwfLGlAYEZEGJyV3K4bhwN8SRitrEGTvcJ3QtN4G477zOnJ2uwgKSu3c98WaM9pM76Rik0i8+X981vdRzi5zUGjAvRvf5tPPx2Ie3uf592sGTNNkzqaDjJm8kGd/2MihonI6Rgfz0Z/O5u6RHb1Wl8KIiDQ4+4q3AtAqoAtGWrLrYEQ7CGzhtZqkKpvVwuRr+xIW4MPafbn8Y87WunkjwyC83y2884cFXOHfGqdh8Er5Pl76bCTli14Dh3fHOjQmG/bnct17y7j9k5XsyiwkMsiXFy9LYtb9QxnZNfqk+8XVB4UREWlwsstdm+N1CO2uLpoGrFV4AK9c2QuAdxbsYuG2zDp7L5/gljx39Uwe7nI9hglfB/tz58a3yP3gfMjeWWfv2xSk55bw0NdrGffGryzZlY2vzcKdIzow75ER3DCwLTar96OA9ysQETlGobELgKTInrB/teugwkiDdGFSLH8c6JpuPeHrtWTm1900XMMwuGng40wZ+ToBFh+WBfhzhTWD3z4eBWs+A61JUsnpNFm15xAv/HcTI16dx7TV+zBN+F3veOY+NJy7z2tNSv4mpm2bxivLX+HPP9/G5uzNXqtXGzyISIOSVZyF05qDaRr0i+0J81a6TiiMNFhPXdydFSmH2Hown1v/tYK/Xt6TpFZhdfZ+I9qex78v/oJH5j1ISsFexkeFcM2vTzJh+ywCx/0TAsLr7L0bMrvDyfKUHH7akM7PG9PJyC8F7Fj80+mUUEDfjiXkO2dwyy87SStMO+75Y/cvpltkt/ovHIUREWlg1hxcB4CzrCVtrKVQcBAMK8T28nJlcjL+Plb+eV1frpi6mHX7crnkn79yed9WPDymC63CA+rkPbu06MLXl05n8qrX+GzLF3wVGsLSvBX89b1z6XXpu9B2cJ28b0NTanfw244sZm1IZ86mgxwqco+hcRAStQbflv+jjMOkAz/trfrcaGsgHYry6FBSRKeycs4+uLueqz9CYUREGpSVFQNWzeI2hGavdR2M6Q6+gd4rSk6rc0wIP90/lFdnb+U/yQf4bs1+flyfxi1DErlrZAdC/X08/p7+Nn8eH/AEwxNG8vTCx9lDDjeGmtz2/R+4I+k2fEZMbLI7PCfvPcwnS3YzZ+NB8o9aCTc80EafzvvZZ0zjYEkqZUCYXxhdIrrQIbwDHUMT6ZixjQ4rPycsv2JT26jOcNFT0O133vkwaKM8EWlgrvvvzazPWYnt0FWs6eGA3yZDvz/BuNe9XZpU07p9h/nLj5tZlpIDQIsgX+47ryPXD2yLTx0NlswtzeWvS15g5p7ZAHQvLWWSpRXtr/zINROrCSh3OJm1IZ0Pf0thTerhyuPRIX5cmBRLp4RsZqe/T3LFhoPhfuHc0esOru5yNb6GFdZ9DfMnweE9rieGtYGRE6HXNWCx1knN1f39rTAiIg2GaZoM+GwwxY4CWuY9xtyQL2H3IvjdP+GsG71dntSAaZr8sjmDST9tZmema1G0xKggHruwK2N6xNTZNNJZKbN48benyXOU4Od08kBeCdedPQFL/z+Bb1CdvGddO1RYxufLU/n3kj2k55UA4GM1GNcrnusGtKFF2GGmJE/hl9RfAPC3+nND9xu4OelmQnyCYfN/Yd5fIHOL6wWDY2DYI65/Uza/Oq1dYUREGp09eXu45LtLMJ02+tin8GnmH6AsH+5cDDE9vF2e1ILd4eTLFXuZ/L9tZBWUAdCvbQT3jOzIiC4t6ySUHCw8yLMLH+W3DNdMrGi7nXPLnAxpM5KBgx8ntJG0lGxNz+ej31L4bs1+SisWlYsK9uX6AW25fmAbnJbDvLvuXaZvn47DdGAxLFze8XLu7H0nMVS0hKz5N2RULKPvHw7nPgjn/Lneuj0VRkSk0ZmxawYTF03EUdSGPwbcxhO7/wQ+QTBxb501I0v9KCi1886Cnby3aBcl5a5frN3iQrlzRAcu7hmH1eLZUGKaJl9t/oLJq/5BobOs8rjVNOnlE8GQDhczpOMldI/sjsVoWKtcbNify8s/beHXHVmAE8PnMO3jCunboZzA4GxS83eTkptCTklO5XNGtB7BA33uoUNWCiR/CltnuXa6Bte/oUF3ux71PNNIYUREGp1Xlr/Cp5s/pSxnMK+Ed+CyPX+BtkPg5pneLk085GBeCe8v2sVny1IpKnMA0C4ykDuGd+CKs1rhZ/Ns6Cyxl7A6fQW/bviU39KWscviqHI+wieEQa2Hcm6rczm31blE+Ed49P1rYu/hLP76v/8xN2U1hl8aFr+D+Phn46TspM/pG92X+ztcSb89q2HtF67ZZ27xfaHvHyHpSgjwzudSGBGRRueGmTeQnJlM8f6r+S5qNz0PfAuD74XRL3m7NPGww0Vl/GvxHj5anMLhiumoMaF+3HZue64b0IYgvzqYBWOaHNjyPb+teIPFBbtZGuBPgeVIq4iBQc+ongxtPZRhrYfRrUW3OulGMk2TzOJMtuRsYVP2JrbkbGFV2noOl2ec8HqbxUa70HYkhiW6HiHtSHQaJObsIXDjf2DvsiMXB0ZCr2uh7/UNomtTYUREGhW7086gzwdR4iihcOcEVkW8R0TuRrjqY+hxubfLkzpSWGrni+WpvL8opXJwZliADzcOasu43vF0ig6um8Gu6esp/3Uy63b8yK8BviwKCGCrn2+VS1r6BHNuzACGdbiIgfGDCfYNrvbLO5wOMooy2Fewj/0F+9mX7/rv/oL97MnbU6WL5WgWRxR9ontwbpvedIroRGJYIq38WmBLWwepS1zBY+8yKMk98iTDCp0ucLWCdBoDNt8TvrY3KIyISKOyNWcrv//v7zGc/pRufZzNAX/GYtrhgfUQ3sbb5UkdK7U7+H7Nft5esIuUrMLK420jA7mgWwznd4+hf9sIz++jcjgVlr8HexZzMGsTi3wNFgUEsCTAn+KjWk1sJiT6hOLjF4rVLxSr1ReLYcFmsWE1rFgsFmyGjVJHKQcKDnCg8AB2p/2kb2sxLIRYWpGd0xJ7SRy28tbcPmAId/ePxa8wzTX9dt9KV/A4kHxk/IebTxC07g8dzoPe10JIrGfvi4cojIhIozJt2zSeW/IcZnFHOu0Zznd+z0JgFDyyA7y4m6jUL4fTZNaGdL5dtZffdmZTVjGLBCAi0IeRXaMZ3T2GoZ1aer4rx2F3TX9NS6Zs/2pWHVzBwuID/Ornw27fmi/aZrPYiA+Kp1VwK1qHtKaVXyStLL7k78tn+9p9hBVnEmfk0CMonw6+h7EVpoHjJONDgmOhzcAjj5iejWJBt+r+/m74n0REmoX1WesBKCtsRW9LxS6srfopiDQzVovBxb3iuLhXHIWldhZtz2T2poPM3ZLBoaJypq/ez/TV+/G1WTi3YxRX90/g/G7RnmkxsdogNglik/Dt+0cGAYMcdh7L2kpqyjz27vsNR/oGnIWZOAxwAA7DwGH1wxHVAWdUF6zhCcQ5DVqXFhNdkIM1/wDOA2txHp6BzV5Y9f3cv4FLKh4AGK51QMJaubZAaDMI2gyA8LZN+t+CwoiINAgbszcC4CxpTX+fRa6D2hyvWQvys3FhUhwXJsVhdzhZuecQczYdZM6mg6TmFDF3SwZzt2TQKjyAPw5sy7VnJxAR5OHxElYbxPSgTUwP2nCPa2fgQ7shZQGkLHQ98jIhbyXsWnnCl7BUPAByzGDSicK3RRvate+MLaI1hLZ2hY/QVhAS16DGfNQXddOIiNeV2EsY+PlAHKaDgu2Ps8j3JVo7D8D106DT+d4uTxoY0zTZdrCA/yTv54vlqZWbw/nZLFzaJ56bBrejR3zd7RpsdzgpKndQXOagqNSO4+BGfFN/JXD/r9gO72avPYwtRaHsNyPZb0aRZkZSHBBH967dGJHUliEdo/D3aR7r5qibRkQajS05W3CYDoKs4WC30tp2wHWi1VlerUsaJsMw6BIbwqMXduW+UZ3479oD/GvJbjbsz+Prlfv4euU+zm4XwU2D2zGmR+wZ7YezO6uQnzemM3vTQXZmFlBU5qgyjuWIrhWPIzpGB3NB9xiu7RZD34RwLB5e2K0pURgREa9zd9FE+nSgo2W362BEIgS28F5R0ij4+1i5qn8Cv+/XmtWph/h48R5+Wp/Git2HWLH7ELGh/ozo0pIe8aF0jw+lW1wogb4n/9Vnmiab0vL4eeNBft6QztaD+Se91moxCPSxEuBrJdDXSqCvjRZBvgzrHMUF3WNJjGqce+F4g8KIiHide/BqgJlIb+Oowasi1WQYBv3atqBf2xYcvLgbny1L5fNlro3lvlyx96jrXBv2dY8LpUd8GD3iQ+kaF8LurCJ+3pjOzxvT2XeouPJ6q8VgYPsWjOkRy4DESIL9bZUBxM9mqbMN/5obhRER8bqNWUcGr55lmec6qDAitRQT6s+ECzpz98gOLNiaybp9uWw8kMvGA3lk5JeyK7OQXZmFzFiXdsLn+/tYGNapJWN6xDKqWzThgc1vQGl9UxgREa/KK8tjd95uAIpyY+ln2eY60WaA94qSJsHPZmV0j1hG9ziyIFhmfikbD+SyKS2PjQfy2HQgj5SsQkL9bZzfLYbRPWIZ3rklAb7NY4BpQ6EwIiJetSnbtb15q+BWBO05SLhRiNMWgCW2l5crk6aoZYgfI7pEM6JLdOWx4jIHPlbD86u7SrUpjIiIV23I2gBAtxbd8SvdCD7giDsLi7XmK16K1IZaQbxPMVBEvModRtoEdaF/RReNrd0gb5YkIvVMYUREvModRiJtHelnbAXASNB4EZHmRGFERLwmqziLg0UHMTAILQki0XLQdSLhbO8WJiL1SmFERLzG3SrSPqw9AQdcX6f5toOACC9WJSL1TWFERLzGHUZ6RPUgJHMVAOnhfbxYkYh4g8KIiHjNhmxXGEmKSiImNxmAgpZa7EykuVEYERGvME2zsmWkZ3hn2pS4ZtI4EwZ6sywR8QKFERHxin0F+8gtzcVmsdG5uBAf7GSaYYTFd/J2aSJSzxRGRMQr3PvRdInogm2fa7zIKmdn4sIDvFmWiHiBwoiIeIW7iyYpKomy3UsAWG12JirYz5tliYgXKIyIiFe4B6/2iOyBbf8KAHYFJGG1aEt2keZGYURE6p3D6ajcIC/JFoKtJIcS04fc8O5erkxEvEFhRETq3a7cXRTbiwmwBdA+Zz8Aa80ORIWFeLkyEfEGhRERqXfu8SLdI7tj3bcccA1ejQ3z92ZZIuIlCiMiUu82Zrtm0iRFJkHqMgBWOjsTpzAi0iwpjIhIvaucSRPSFrK3A7Da2YmYUIURkeZIYURE6lWZo4yth7YC0KO0FIA9RmsOE0JcmNYYEWmOFEZEpF5tzdmK3Wkn3C+c1gddS8Avd7hWXY1Vy4hIs6QwIiL1qnJ9kageGBWDV91hJDpUC56JNEcKIyJSryrHi0R0g/2rAVjp7EJkkC/+PlZvliYiXqIwIiL1yr0nTZIlAByllPm1IMWM1eBVkWZMYURE6k1heSG7cncBkJSXA8DBsN6AoWm9Is2YwoiI1JtN2ZswMYkJjCEqbT0AKQE9AYhRGBFpthRGRKTeuMeL9IxKgr2uxc7WW7oAEKduGpFmS2FEROqNO4z0CIyHwkyw+rLa3g5AS8GLNGMKIyJSb9ZnubpmksrsrgPxfdmb5wAURkSaM4UREakX6YXppBWmYTEs9MzZ5zqYMIC03BIADWAVacYURkSkXqzNXAtA54jOBO1bBUBJ3Nnkl7haSWK1FLxIs6UwIiL1wh1Gekd0hcwtAKSF9QYg2M9GsJ/Na7WJiHcpjIhIvVib4QojfYxA14HIjqSVBQEaLyLS3CmMiEidK7GXsClnEwC9Cw65Dh41XkQb5Ik0bwojIlLnNmVvwu60E+kfSesDrum9JAwgPa8ijKhlRKRZUxgRkTqXnJkMQJ+WvTAqNsejzUDSNZNGRFAYEZF6UDlexC8a7MUQEAGRnSq7abRJnkjzpjAiInXKNM0jLSOlZa6DCQPAYuFgnlpGRERhRETq2L78feSU5GCz2OiW6dqxl4QBAGoZERGglmFk6tSpJCYm4u/vT79+/Vi0aNFJr50/fz6GYRz32LJlS62LFpHGw90q0r1Fd/x2/+o62HYwZXYnWQWlgFpGRJq7GoeRr776igceeIAnn3ySNWvWMHToUMaOHUtqauopn7d161bS0tIqH506dap10SLSeCRnJAPQJyAGig+Bfzi06k9GvqtVxNdqoUWQr/cKFBGvq3EYee2117j11lu57bbb6NatG5MnTyYhIYG33nrrlM+Ljo4mNja28mG1WmtdtIg0Hu6VV/sUFroOdBwFVlvlTJqYMD8Mw/BWeSLSANQojJSVlbFq1SpGjx5d5fjo0aNZvHjxKZ/bt29f4uLiGDVqFPPmzTvltaWlpeTl5VV5iEjjU1BWwPbD2wHovd+1Yy+dXD8/3GuMxIVqTxqR5q5GYSQrKwuHw0FMTEyV4zExMaSnp5/wOXFxcbz77rtMmzaN6dOn06VLF0aNGsXChQtP+j6TJk0iLCys8pGQkFCTMkWkgViftR6n6SQ+MIbo9I2AAR3PBziqZUTjRUSau1rtTHVsk6ppmidtZu3SpQtdunSp/H7QoEHs3buXV199lWHDhp3wORMnTmTChAmV3+fl5SmQiDRC7sGrvX0iXAda9YOgKODITBoNXhWRGrWMREVFYbVaj2sFycjIOK615FQGDhzI9u3bT3rez8+P0NDQKg8RaXwqFzsrrOhq7Tym8lzlUvCa1ivS7NUojPj6+tKvXz/mzJlT5ficOXMYPHhwtV9nzZo1xMXF1eStRaSRcZpO1mWuA6DPgc2ug52OjDdzd9NoXxoRqXE3zYQJE7jhhhvo378/gwYN4t133yU1NZXx48cDri6W/fv388knnwAwefJk2rVrR48ePSgrK+PTTz9l2rRpTJs2zbOfREQalF2Hd5Ffnk+AxZfOhbkQHAOxvSrPK4yIiFuNw8g111xDdnY2L7zwAmlpaSQlJTFz5kzatm0LQFpaWpU1R8rKynj44YfZv38/AQEB9OjRgx9//JGLLrrIc59CRBoc93iRJGuI6wdNpwvA4mqMdTrNyqXg1U0jIoZpmqa3izidvLw8wsLCyM3N1fgRkUbiqV+f4j87/8PtpRbuO7Abrv4Eul8KQEZ+Cef85RcsBmx9aSw+Vu1MIdIUVff3t34CiEidqFzs7FA6WGzQfmTluYO5rmXgW4b4KYiIiMKIiHjeoZJD7M7bDUCv0jJoMwj8j/xfUVpuMaAuGhFxURgREY9zz6JJxIdwp7PKlF7gyHgRDV4VERRGRKQOVC52VpDrOtCpahg5suCZloIXEYUREakDlTv1FhdDeFuIqrpLd+VS8OqmEREURkTEw8qd5WzI2gBAn9JSVxfNMdtFVG6Sp24aEUFhREQ8bNuhbZQ4SghxmiSW26usuuqmlhEROZrCiIh4lLuLpndJCRZbALQ7t8p50zTVMiIiVSiMiIhHuTfH611aConDwKfqINW8EjtFZQ5As2lExEVhREQ8yj2Tpk9JKXQ+eRdNeKAP/j7W+ixNRBoohRER8ZiDhQdJK0zDYpr0LC078XgR7UkjIsdQGBERj3EvAd+5rJygqK4Q3ua4a9Ldq6+qi0ZEKiiMiIjHVC52VnriLhqA9Ip9aTR4VUTcFEZExGPWuhc7Kyk9YRcNQHqee18arb4qIi4KIyLiESX2EjZlbwKgt+kLCQNOeJ17KfjYML96q01EGjaFERHxiE3Zm7CbDiLtDlq3HQ5WnxNel14ZRtQyIiIuCiMi4hHuwat9Sksxulx40us0m0ZEjqUwIiIekXxgGVAxXqTjBSe8pqTcweGickCzaUTkCIURETljpmmSnLkGgD4h7SC45Qmvc3fRBPpaCfW31Vd5ItLAKYyIyBnblLOJHHsRAU4n3TqOPel1lYNXQ/0xjtnJV0SaL4URETljc1NmAzCkuAS/zicfL3LQPV5EXTQichSFERE5Y/N2/QjAeWYAxPU96XVHt4yIiLgpjIjIGdmbv5ftxQexmibDOl0KlpP/WFHLiIiciMKIiJyReTtcrSL9SkoJ63PDKa9Nq9iXRkvBi8jRFEZE5IzM3fEdAOfZIiC66ymv3Zvj3iRPC56JyBEKIyJSa4dKDrGm8AAAIztfccprC0vtbD2YD0BSq9A6r01EGg+FERGptQVbvsFpQNfSMuL73nzKa1enHsLhNGkdEUCcWkZE5CgKIyJSa/O2ubpoRvrFQkjMKa9dnpIDwDntWtR5XSLSuCiMiEitFJcXsbhoHwDndbnytNe7w8jZiQojIlKVwoiI1MqS9f+mxIB4u4MufW855bWldgfJew8DcI7CiIgcQ2FERGpl3rbpAIwMaI3hF3zKazfsz6XU7iQq2Jf2UUH1UZ6INCIKIyJSY46yYhYUu7poRna+/LTXL6voounftoX2pBGR4yiMiEiNJSd/yCGLhVCnyVl9Tj2LBmCFxouIyCkojIhIjc3dPg2A4YEJ+NhOvZqqw2mycs8hQDNpROTEFEZEpEbMokPMK04DYGQ1ZtFsTc8nv8ROsJ+NbnEhdV2eiDRCCiMiUiM71nzAXh8bviYM6f6H016/Yreri+asthHYrPqRIyLH008GEamRedu+B2BgUAKBvqefGXNksbOIuixLRBoxhRERqb7DqcwtzwDgvNPsRQNgmibLK1pGztZ4ERE5CYUREam29NUfsdHPD8OE4Z0vO+31e7KLyMwvxddqoXdCeJ3XJyKNk8KIiFSPaTJ/+/cA9A5qRVRA1Gmf4m4V6dU6DH8fa11WJyKNmMKIiFRPWjJzzXwARna8tFpPqRwvovVFROQUFEZEpFrykj9jhb9rTZHz2o+t1nPcM2m02JmInIrCiIicnsPOrzv+i90waB8QTbuwdqd9SkZeCXuyizAM6NdWM2lE5OQURkTk9HbNZ561HICR7S+u1lPc40W6xYYS6u9TZ6WJSOOnMCIip1W29nMWBQYAcF7b86v1nBUaLyIi1aQwIiKnVpLHit1zKLRYaOkbTlJUUrWetkxhRESqSWFERE7t19eY6+ealjui3QVYjNP/2MgtLmfrQdfMGy12JiKnozAiIieXvRPnkjeZX9FFMzJhZLWetmpPDqYJiVFBtAzxq8sKRaQJUBgRkZObNZEFflYybDYCbYEMiBtQractTzkEwNnaj0ZEqkFhRERObNvPlOyYzSstXIHimq7X4Gv1rdZTl6dkA+qiEZHqURgRkePZS2HWRD4MC2W/j43owGjG9xpfraeWlDtYvz8XgAGJkXVZpYg0EQojInK8pVPZm7eHD8JDAXj07EcJ9Ams1lPXpB6m3GESE+pHQouAuqxSRJoIhRERqSovDXPB35kUGUGZYTAobhCj246u9tMrl4Bv1wLDMOqqShFpQhRGRKSq/z3LPB8niwIDsFlsTBwwsUahwh1GtL6IiFSXwoiIHJG6lOL1X1cOWv1Tjz+RGJZY7afbHU5W7XHPpFEYEZHqURgRERenA2Y+wvthoRzwsREXFMftPW+v0UtsPJBHUZmDUH8bXWJC6qhQEWlqFEZExGX1J+zJ2sRHFYNWHzv7sWoPWnU7eryIxaLxIiJSPQojIgJFOZi/vMBfIyMoNwyGtBrCeW3Oq/HLLK/Yj+ZsjRcRkRpQGBERmD+JX4xiFgcG4GPxYeI5NRu0CuB0mlVaRkREqkthRKS5O7iRopUf8Eqka9DqzUk30za0bY1fZmdmAYeKyvH3sdCzVZinqxSRJkxhRKQ5M02Y+SjvhQWTbrMRHxTPbT1vq9VLLa9oFemTEI6vTT9aRKT69BNDpLnKT4fpf2bXgaV8HOYatPr4OY8TYKvdqqkrUtzri2gJeBGpGZu3CxCReuYoh2XvwPyXMcvymRQbjd0wGNZ6GCMSRtTqJQ/mlfDbTtfmeOdovIiI1JDCiEhzkrIQZj4CmVtIt1p5rU0nllpL8bX48vjZj9dq+fbFO7K478s1ZBWUERXsx1ltwz1ft4g0aQojIs1B7n6Y/SRs/I4iw+CDlrH8KySQUrMUA4MH+z1IQmhCjV7S6TR5c94O/u9/23Ca0DU2hKnXn0Wgr36siEjN6KeGSFNmL4Ulb8LCv+MsL+KHkGCmRMeR6SwF006/mH48evajdI/sXqOXzSks48GvklmwLROAq/u35oVLk/D3sdbFpxCRJq5WA1inTp1KYmIi/v7+9OvXj0WLFp3y+gULFtCvXz/8/f1p3749b7/9dq2KFZFqME3I2AyL34C3BsMvz7PS6uDadh14OqoFmc5SWge35v9G/B8fjfmoxkFkdeohLpmyiAXbMvGzWfjb73vxt9/3VhARkVqrccvIV199xQMPPMDUqVMZMmQI77zzDmPHjmXTpk20adPmuOtTUlK46KKLuP322/n000/57bffuOuuu2jZsiVXXnmlRz6ESLNXmAW75sPOua5HfhoAe202Xotrxf/8rUA5wT7B3NHrDq7rdh2+Vt8avYVpmnz0227+OnMzdqdJYlQQU68/i25xoZ7/PCLSrBimaZo1ecKAAQM466yzeOuttyqPdevWjcsuu4xJkyYdd/1jjz3GDz/8wObNmyuPjR8/nrVr17JkyZJqvWdeXh5hYWHk5uYSGqoffNKMOZ1QVuB6ZO+EnXMp3/kL6Zkb2Wuzsc/Hxj6bjf2+vuwLDGMb5dhxYjEsXNX5Ku7qcxct/Gs+2yW/pJxHv13HTxvSAbi4ZxwvX9mTEH8fT39CEWlCqvv7u0YtI2VlZaxatYrHH3+8yvHRo0ezePHiEz5nyZIljB49usqxMWPG8MEHH1BeXo6Pz/E/zEpLSyktLa3yYerC5K/vY3fe5tNfKOIxJpgATtfXmJiYGGbFf3FWXObAMO0Y2MG0u77HjokTE3AAhRYL+2020vysOBPiT/Bern9DMbZe9Aq8AUdmAv+ckwakVa3IhFK7k9JyB8XlDkrKHZSUOymxOyguc1Bqd5JdUEpeiR0fq8GTF3XjpsHtajXzRkTkRGoURrKysnA4HMTExFQ5HhMTQ3p6+gmfk56efsLr7XY7WVlZxMXFHfecSZMm8fzzz9ektFrZcGgFy/wL6vx9RM6Mgeuf6in+uTptOMpbYJa1wFneAmfFf82yKHaUtWQHDmD3GVXRKjyAN67rS982EWf0OiIix6rVbJpj/4/INM1T/l/Sia4/0XG3iRMnMmHChMrv8/LySEio2bTD6jir5VBaHN7k8dcVOTUDAwsYBiYWXGHDgIpjYGAaNpwWP5yGD6bFH4fFF6fhh9Pih8Pij2n4YrX4EWKNIdgaQ4AlDMOo/YLKfjYr/j4W/H2sRx62I98H+FjpFBOsQaoiUidqFEaioqKwWq3HtYJkZGQc1/rhFhsbe8LrbTYbkZEnXjbaz88PPz+/mpRWK3dd8bc6fw8RERE5tRr9r5Svry/9+vVjzpw5VY7PmTOHwYMHn/A5gwYNOu762bNn079//xOOFxEREZHmpcbtuhMmTOD999/nww8/ZPPmzTz44IOkpqYyfvx4wNXFcuONN1ZeP378ePbs2cOECRPYvHkzH374IR988AEPP/yw5z6FiIiINFo1HjNyzTXXkJ2dzQsvvEBaWhpJSUnMnDmTtm3bApCWlkZqamrl9YmJicycOZMHH3yQN998k/j4eKZMmaI1RkRERASoxToj3qB1RkRERBqf6v7+rv3wexEREREPUBgRERERr1IYEREREa9SGBERERGvUhgRERERr1IYEREREa9SGBERERGvUhgRERERr1IYEREREa+q8XLw3uBeJDYvL8/LlYiIiEh1uX9vn26x90YRRvLz8wFISEjwciUiIiJSU/n5+YSFhZ30fKPYm8bpdHLgwAFCQkIwDMNjr5uXl0dCQgJ79+7Vnjf1QPe7ful+1y/d7/ql+13/anPPTdMkPz+f+Ph4LJaTjwxpFC0jFouF1q1b19nrh4aG6i9zPdL9rl+63/VL97t+6X7Xv5re81O1iLhpAKuIiIh4lcKIiIiIeFWzDiN+fn48++yz+Pn5ebuUZkH3u37pftcv3e/6pftd/+rynjeKAawiIiLSdDXrlhERERHxPoURERER8SqFEREREfEqhRERERHxqmYdRqZOnUpiYiL+/v7069ePRYsWebukJmHhwoWMGzeO+Ph4DMPg+++/r3LeNE2ee+454uPjCQgIYMSIEWzcuNE7xTYBkyZN4uyzzyYkJITo6Gguu+wytm7dWuUa3XPPeeutt+jVq1flwk+DBg3ip59+qjyve113Jk2ahGEYPPDAA5XHdL8967nnnsMwjCqP2NjYyvN1db+bbRj56quveOCBB3jyySdZs2YNQ4cOZezYsaSmpnq7tEavsLCQ3r1788Ybb5zw/N/+9jdee+013njjDVasWEFsbCwXXHBB5R5EUjMLFizg7rvvZunSpcyZMwe73c7o0aMpLCysvEb33HNat27Nyy+/zMqVK1m5ciXnnXcel156aeUPZN3rurFixQreffddevXqVeW47rfn9ejRg7S0tMrH+vXrK8/V2f02m6lzzjnHHD9+fJVjXbt2NR9//HEvVdQ0AeZ3331X+b3T6TRjY2PNl19+ufJYSUmJGRYWZr799tteqLDpycjIMAFzwYIFpmnqnteHiIgI8/3339e9riP5+flmp06dzDlz5pjDhw8377//ftM09Xe7Ljz77LNm7969T3iuLu93s2wZKSsrY9WqVYwePbrK8dGjR7N48WIvVdU8pKSkkJ6eXuXe+/n5MXz4cN17D8nNzQWgRYsWgO55XXI4HHz55ZcUFhYyaNAg3es6cvfdd3PxxRdz/vnnVzmu+103tm/fTnx8PImJiVx77bXs2rULqNv73Sg2yvO0rKwsHA4HMTExVY7HxMSQnp7upaqaB/f9PdG937NnjzdKalJM02TChAmce+65JCUlAbrndWH9+vUMGjSIkpISgoOD+e677+jevXvlD2Tda8/58ssvWb16NStWrDjunP5ue96AAQP45JNP6Ny5MwcPHuSll15i8ODBbNy4sU7vd7MMI26GYVT53jTN445J3dC9rxv33HMP69at49dffz3unO6553Tp0oXk5GQOHz7MtGnTuOmmm1iwYEHled1rz9i7dy/3338/s2fPxt/f/6TX6X57ztixYyu/7tmzJ4MGDaJDhw7861//YuDAgUDd3O9m2U0TFRWF1Wo9rhUkIyPjuMQnnuUela1773n33nsvP/zwA/PmzaN169aVx3XPPc/X15eOHTvSv39/Jk2aRO/evXn99dd1rz1s1apVZGRk0K9fP2w2GzabjQULFjBlyhRsNlvlPdX9rjtBQUH07NmT7du31+nf72YZRnx9fenXrx9z5sypcnzOnDkMHjzYS1U1D4mJicTGxla592VlZSxYsED3vpZM0+See+5h+vTpzJ07l8TExCrndc/rnmmalJaW6l572KhRo1i/fj3JycmVj/79+3P99deTnJxM+/btdb/rWGlpKZs3byYuLq5u/36f0fDXRuzLL780fXx8zA8++MDctGmT+cADD5hBQUHm7t27vV1ao5efn2+uWbPGXLNmjQmYr732mrlmzRpzz549pmma5ssvv2yGhYWZ06dPN9evX2/+4Q9/MOPi4sy8vDwvV9443XnnnWZYWJg5f/58My0trfJRVFRUeY3uuedMnDjRXLhwoZmSkmKuW7fOfOKJJ0yLxWLOnj3bNE3d67p29Gwa09T99rSHHnrInD9/vrlr1y5z6dKl5iWXXGKGhIRU/m6sq/vdbMOIaZrmm2++abZt29b09fU1zzrrrMqpkHJm5s2bZwLHPW666SbTNF3Tw5599lkzNjbW9PPzM4cNG2auX7/eu0U3Yie614D50UcfVV6je+45t9xyS+XPjZYtW5qjRo2qDCKmqXtd144NI7rfnnXNNdeYcXFxpo+PjxkfH29eccUV5saNGyvP19X9NkzTNM+sbUVERESk9prlmBERERFpOBRGRERExKsURkRERMSrFEZERETEqxRGRERExKsURkRERMSrFEZERETEqxRGRERExKsURkRERMSrFEZERETEqxRGRERExKsURkRERMSr/h+sXdk7pYn+hQAAAABJRU5ErkJggg==", "text/plain": [ "
" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "plot.plot(classical_rdf, label=\"Classical\")\n", "plot.plot(rpmd_rdf, label=\"RPMD\")\n", "plot.plot(qtb_rdf, label=\"adQTB\")\n", "plot.legend()\n", "plot.show()" ] }, { "cell_type": "markdown", "id": "00431d18-a236-4b82-8fc3-6b1234468e8d", "metadata": {}, "source": [ "The results from RPMD and adQTB are very similar, while the classical simulation produced a quite different distribution. That shows the importance of NQEs in this system.\n", "\n", "## When Should You Include Nuclear Quantum Effects?\n", "\n", "You might be tempted to add nuclear quantum effects to all your simulations, especially given the relatively low cost of adQTB, but that would not be correct. Most force fields are parameterized at least partly based on experimental data. That means they already implicitly take NQEs into account. They are designed to match experiment *without* explicitly adding NQEs.\n", "\n", "On the other hand, quantum chemistry methods like Density Functional Theory treat nuclei classically, using quantum mechanics only for the electrons. QM/MM simulations can therefore benefit from adding explicit NQEs. The same is true of force fields and machine learning potentials that are parameterized solely from quantum chemistry data. The potential function reflects a purely classical model of the nuclei. Using RPMD or adQTB to add quantum effects increases the realism of the model." ] } ], "metadata": { "kernelspec": { "display_name": "Python 3 (ipykernel)", "language": "python", "name": "python3" }, "language_info": { "codemirror_mode": { "name": "ipython", "version": 3 }, "file_extension": ".py", "mimetype": "text/x-python", "name": "python", "nbconvert_exporter": "python", "pygments_lexer": "ipython3", "version": "3.13.3" } }, "nbformat": 4, "nbformat_minor": 5 }