{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "# Quickstart\n", "\n", "A quick introduction on how to use the OQuPy package to compute the dynamics of a quantum system that is possibly strongly coupled to a structured environment. We illustrate this by applying the TEMPO method to the strongly coupled spin boson model.\n", "\n", "- [launch binder](https://mybinder.org/v2/gh/tempoCollaboration/OQuPy/HEAD?labpath=tutorials%2Fquickstart.ipynb) (runs in browser),\n", "- [download the jupyter file](https://raw.githubusercontent.com/tempoCollaboration/OQuPy/main/tutorials/quickstart.ipynb), or\n", "- read through the text below and code along." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "**Contents:**\n", "\n", "* Example - The spin boson model\n", " * 1. The model and its parameters\n", " * 2. Create system, correlations and bath objects\n", " * 3. TEMPO computation" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "First, let's import OQuPy and some other packages we are going to use" ] }, { "cell_type": "code", "execution_count": 1, "metadata": {}, "outputs": [], "source": [ "import sys\n", "sys.path.insert(0,'..')\n", "\n", "import oqupy\n", "import numpy as np\n", "import matplotlib.pyplot as plt" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "and check what version of tempo we are using." ] }, { "cell_type": "code", "execution_count": 2, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "'0.4.0'" ] }, "execution_count": 2, "metadata": {}, "output_type": "execute_result" } ], "source": [ "oqupy.__version__" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Let's also import some shorthands for the spin Pauli operators and density matrices." ] }, { "cell_type": "code", "execution_count": 3, "metadata": {}, "outputs": [], "source": [ "sigma_x = oqupy.operators.sigma(\"x\")\n", "sigma_y = oqupy.operators.sigma(\"y\")\n", "sigma_z = oqupy.operators.sigma(\"z\")\n", "up_density_matrix = oqupy.operators.spin_dm(\"z+\")\n", "down_density_matrix = oqupy.operators.spin_dm(\"z-\")" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "-------------------------------------------------\n", "## Example - The spin boson model\n", "As a first example let's try to reconstruct one of the lines in figure 2a of [Strathearn2018] ([Nat. Comm. 9, 3322 (2018)](https://doi.org/10.1038/s41467-018-05617-3) / [arXiv:1711.09641v3](https://arxiv.org/abs/1711.09641)). In this example we compute the time evolution of a spin which is strongly coupled to an ohmic bath (spin-boson model). Before we go through this step by step below, let's have a brief look at the script that will do the job - just to have an idea where we are going:" ] }, { "cell_type": "code", "execution_count": 4, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "--> TEMPO computation:\n", "100.0% 150 of 150 [########################################] 00:00:03\n", "Elapsed time: 3.9s\n" ] }, { "data": { "text/plain": [ "" ] }, "execution_count": 4, "metadata": {}, "output_type": "execute_result" }, { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAY4AAAEICAYAAABI7RO5AAAAOXRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjMuNCwgaHR0cHM6Ly9tYXRwbG90bGliLm9yZy8QVMy6AAAACXBIWXMAAAsTAAALEwEAmpwYAAAoCklEQVR4nO3deXhV5b328e8ve2cgA4QMjAkmQBiVQQJWBepctVZqnbA9byetnlZr+3bSnp637WnP6bGtHY92sHZunU8HWhUnahHrQFRQmSTMCWACIUACmX/vH9nYCEGzIclae+f+XBdX9l57kdwJCXfW86z1LHN3REREeiol6AAiIpJYVBwiIhIXFYeIiMRFxSEiInFRcYiISFyiQQfoawUFBV5SUhJ0DBGRhPLCCy/scvfC7l5L+uIoKSmhoqIi6BgiIgnFzLYc7TUNVYmISFxUHCIiEhcVh4iIxCXp5zhERFpbW6mqqqKpqSnoKKGTkZFBUVERqampPf47Kg4RSXpVVVXk5ORQUlKCmQUdJzTcnd27d1NVVUVpaWmP/16ohqrM7HwzW2dmlWZ2czevf9jMas1sRezPNUHkFJHE0tTURH5+vkrjMGZGfn5+3EdioTniMLMIcDtwLlAFLDezRe6++rBd73X3G/o9oIgkNJVG947l6xKmI445QKW7b3T3FuAeYEFQYeoaW/ivB1ez/vX9QUUQEQmlMBXHaGBbl+dVsW2Hu9TMXjazB8ysuLt3ZGbXmlmFmVXU1tYeU5gXtuzhl09v5tzvLeXSH/+DbXUHjun9iIgkmzAVR0/8BShx92nAY8Cvu9vJ3e9w93J3Ly8s7PaK+bd17pThPPtvZ/OlCyez/vX9fPAXz1PX2HLsyUVEkkSYiqMa6HoEURTb9gZ33+3uzbGndwKz+jJQQXY6H5s/lp9/eDbb6w/y0V8t52BLe19+SBGRt7V48WImTpzI+PHjueWWW7rdp6mpiTlz5jB9+nSmTp3KV77ylV77+GEqjuVAmZmVmlkasBBY1HUHMxvZ5enFwJr+CDa7JI8fLJzJim31/PjvG/rjQ4qIdKu9vZ3rr7+ehx9+mNWrV3P33XezevXh5xBBeno6S5YsYeXKlaxYsYLFixfz7LPP9kqG0BSHu7cBNwCP0FkI97n7KjP7mpldHNvtRjNbZWYrgRuBD/dXvvNPHMF7po/ip3/fQHX9wf76sCKSRFatWsU555zDhAkT+PrXv84nP/lJli9fHtf7eP755xk/fjxjx44lLS2NhQsX8uc///mI/cyM7OxsoPMCyNbW1l47syw0p+MCuPtDwEOHbftyl8dfBL7Y37kOufmCSTy6aifffHgtP7xqZlAxROQ4/MdfVrF6+75efZ9TRg3mK++Z+pb7NDU1cfnll3P//fczduxYJk2axKxZs5g9e/Yb+8ybN4/9+488k/PWW2/lnHPOAaC6upri4n+O6hcVFfHcc891+zHb29uZNWsWlZWVXH/99ZxyyinH8ukdIVTFEXajcwdx3fyx/HBJJR+dW8qM4tygI4lIgnj88ceZOXMmU6d2FkxLSwuf/exn37TPU0891asfMxKJsGLFCurr67nkkkt49dVXOfHEE4/7/ao44nTdO8fxq39s5mdLN3L7B04OOo6IxOntjgz6yooVK5g5s3OkYvv27WRnZ3P66ae/aZ+eHHGMHj2abdv+eeVCVVUVo0d3d+XCP+Xm5nLmmWeyePFiFUcQstKjXHXKGH62dCPb6g5QnJcZdCQRSQBpaWlUV3eeKPrFL36RlpYjT+/vyRHH7NmzWb9+PZs2bWL06NHcc8893HXXXUfsV1tbS2pqKrm5uRw8eJDHHnuMm2666fg/EUI0OZ5IPnRq50Jpv/7H5qCjiEiCeP/738/SpUuZOHEi06dP59RTT+XTn/503O8nGo1y22238a53vYvJkydzxRVXvDH8BXDhhReyfft2duzYwZlnnsm0adOYPXs25557LhdddFGvfC7m7r3yjsKqvLzc++LWsZ+8+yWeXFvDM/92NtnpOnATCbM1a9YwefLkoGOEVndfHzN7wd3Lu9tfRxzH6Oq5pexvbuOPL1W//c4iIklExXGMphcNYdKIHB54oSroKCIi/UrFcYzMjMtmFbFyWz2VNVpBVyTskn1Y/lgdy9dFxXEcFswYTSTFeOAFDVeJhFlGRga7d+9WeRzm0B0AMzIy4vp7mtU9DoU56ZwxoZA/vlTF5981kUiKbhQjEkZFRUVUVVVxrLdZSGaH7jkeDxXHcbpsVhFPrK1hWeUu3jnh2JZwF5G+lZqaGtc9teWtaajqOJ01eRg56VH+unJ70FFERPqFiuM4pUcjnDtlOI+s2klLW0fQcURE+pyKoxe8e9pI9jW18fSGXUFHERHpcyqOXjC3rICc9CgPvrwj6CgiIn1OxdEL0qMRzp06nEc1XCUiA4CKo5e8+6TYcFWlhqtEJLmpOHrJ3LICstIiPLr69aCjiIj0KRVHL0mPRphXVsiSta/r6lQRSWoqjl509uRhvL6vmVere/d+xiIiYaLi6EVnThqGGTy+RsNVIpK8VBy9qCA7nZnFuSxZWxN0FBGRPqPi6GVnTx7OK9V7eX1fU9BRRET6hIqjl50zeTgAT6zRUYeIJKdQFYeZnW9m68ys0sxufov9LjUzN7Nu74cbpAnDsykaOognNM8hIkkqNMVhZhHgduACYApwlZlN6Wa/HOBTwHP9m7BnzIxzJg9nWeUuDra0Bx1HRKTXhaY4gDlApbtvdPcW4B5gQTf7fR34JhDaSYSzJw+jua1DV5GLSFIKU3GMBrZ1eV4V2/YGMzsZKHb3B9/qHZnZtWZWYWYVQdzx65TSfLLTozyxVsNVIpJ8wlQcb8nMUoDvAp99u33d/Q53L3f38sLC/r8rX1o0hfkTCnhiTY2uIheRpBOm4qgGirs8L4ptOyQHOBF40sw2A+8AFoVxghzgrEnDqdmvq8hFJPmEqTiWA2VmVmpmacBCYNGhF919r7sXuHuJu5cAzwIXu3tFMHHf2pkTCzFDw1UiknRCUxzu3gbcADwCrAHuc/dVZvY1M7s42HTxy89OZ9roITy1XhPkIpJcokEH6MrdHwIeOmzbl4+y7xn9kel4zJ9QyI+e3MDeg60MGZQadBwRkV4RmiOOZDR/QiHtHc4/dFquiCQRFUcfmlGcS056lKXr+/+UYBGRvqLi6EOpkRROG5/P0td26bRcEUkaKo4+Nn9CIdX1B9lQ2xh0FBGRXqHi6GPzyzovQFz6moarRCQ5qDj6WHFeJqUFWZrnEJGkoeLoB/PLCnh2426aWrVarogkPhVHP5g/oZCm1g4qNu8JOoqIyHFTcfSDd4zNJzViGq4SkaSg4ugHWelRyk/I0wS5iCQFFUc/mT+hkLU79/P6vtDef0pEpEdUHP1k/oQCQKflikjiU3H0kykjB1OQna7VckUk4ak4+omZMXd8Pk9X7qKjQ8uPiEjiUnH0o7llhexubGHNTt0VUEQSl4qjH80r65zn0HCViCQyFUc/Gj44gwnDs1mm4hCRBKbi6Gdzxxfy/OY6LT8iIglLxdHP5k0ooKWtg+c31QUdRUTkmKg4+tkppXmkRVJYptvJikiCUnH0s8y0KCefkKsJchFJWCqOAMwrK2TNjn3U7m8OOoqISNxUHAE4dFru0xquEpEEFKriMLPzzWydmVWa2c3dvP6vZvaKma0ws2VmNiWInMdr6qgh5GamarhKRBJSaIrDzCLA7cAFwBTgqm6K4S53P8ndZwDfAr7bvyl7RyTFOH18AU+tr8Vdy4+ISGIJTXEAc4BKd9/o7i3APcCCrju4e9e1OrKAhP1fd974Amr2N7O+piHoKCIicQlTcYwGtnV5XhXb9iZmdr2ZbaDziOPG7t6RmV1rZhVmVlFbG85lzOdq+RERSVBhKo4ecffb3X0ccBPw70fZ5w53L3f38sLCwv4N2ENFQzMZW5DFU7qdrIgkmDAVRzVQ3OV5UWzb0dwDvLcvA/W1uWUFPLexjuY2LT8iIokjTMWxHCgzs1IzSwMWAou67mBmZV2evhtY34/5et28skIOtrbz4pb6oKOIiPRYaIrD3duAG4BHgDXAfe6+ysy+ZmYXx3a7wcxWmdkK4DPAh4JJ2zveMTaPSIppuEpEEko06ABduftDwEOHbftyl8ef6vdQfSgnI5WZxbksq9zFF4IOIyLSQ6E54hio5pUV8kr1XvY0tgQdRUSkR1QcAZtbVoA7PL1Bp+WKSGJQcQRsetEQcjKiuiugiCQMFUfAopEUThuXz1Prd2n5ERFJCCqOEJhbVkh1/UE27WoMOoqIyNtScYTA/NjyI7oroIgkAhVHCJyQn0Vx3iCtWyUiCUHFERJzxxfyzIbdtLZ3BB1FROQtqThCYn5ZAQ3NbazcVh90FBGRt6TiCInTxhWQYlpmXUTCT8UREkMyUzmpKFfrVolI6Kk4QmR+WQErq/ayr6k16CgiIkel4giRueMLaO9wntmwO+goIiJHpeIIkZljhpKZFtFwlYiEmoojRNKiKZw6Nl/rVolIqKk4QmZuWQGbdx9gW92BoKOIiHRLxREy82LLj+i0XBEJKxVHyIwrzGbkkAyWVWqeQ0TCScURMmbG3PEFPF25m/YOLbMuIuGj4giheRMK2XuwlVeq9wYdRUTkCCqOEDp9XD4Ay3RaroiEkIojhPKz05k6arAmyEUklHpUHGY21Mw+aGZ/NLNVZvZXM/uYmQ3r64AD1byyQl7cuofG5rago4iIvMnbFoeZ/QH4AzAMuMndpwKfADKA35nZk32acICaV1ZAa7vz3CYtPyIi4dKTI46PuvuZ7n6ru78G4O5b3f1/3P084L29FcbMzjezdWZWaWY3d/P6Z8xstZm9bGZPmNkJvfWxw2bWCUNJj6ZouEpEQudti8Pd6w/fZma/NbOM3gxiZhHgduACYApwlZlNOWy3l4Byd58GPAB8qzczhElGaoRTxuarOEQkdI51ctyAn5jZIOAzvZRlDlDp7hvdvQW4B1jQdQd3/5u7H1qL41mgqJc+dijNG19AZU0DO/YeDDqKiMgbjrU4NgFfBX4MZPVSltHAti7Pq2LbjuZq4OHuXjCza82swswqamsT95TWubHlR7TooYiEybEWxx3uvpnO8ji/19L0kJn9C1AOfLu71939Dncvd/fywsLC/g3XiyaNyKEgO13DVSISKj0uDjP76qHH7r4t9nZz7Cyr3lANFHd5XhTbdniOc4AvARe7e3MvfexQMjPmlRXwdOUuOrT8iIiERDxHHF82s2+a2c/M7ONmNrSXsywHysys1MzSgIXAoq47mNlM4Kd0lkZNL3/8UJo7voDdjS2s2bkv6CgiIkB8xeFAE/AInUcG/zCz6b0VxN3bgBti738NcJ+7rzKzr5nZxbHdvg1kA/eb2QozW3SUd5c0tMy6iISNufdsCMTMVnUdljKzCcBP3P2svgrXG8rLy72ioiLoGMfl/O8vJS8rjbs+9o6go4jIAGFmL7h7eXevxXPEscvMZh16ErsYMHFnnhPIOycWsnxzHQ1afkREQiCe4riRziVGfmdmN5nZXXSelit97KyJw2htd62WKyKh0JO1qlIA3H0lMAO4O/bSEmDhodel78w6YSg5GVGWrB0Q5wOISMj15D/9x8zsXjO7Ckh39weB/wH20HmG04t9GVAgGklh/oRC/rauVqflikjgerJW1dnAfwAlwINm9izwBHAS8D13n9GXAaXTWROHUbu/mVXbdVquiAQr2pOd3H01sBr4bzMb5O5aPKmfnTGxEDNYsraGk4qGBB1HRAawuOcnVBrByM9OZ3pRLkvWaZ5DRIKlie0EctakYbxcVc+uhqReaUVEQk7FkUDOmjQMd/j7Op2WKyLBUXEkkKmjBjMsJ13DVSISKBVHAjEzzpw4jKWv1dLa3hF0HBEZoFQcCebMScPY39TGC1v2BB1FRAYoFUeCmVtWQGrE+JuuIheRgKg4Ekx2epRTSvN5fM3rQUcRkQFKxZGAzps6nA21jWyobQg6iogMQCqOBHTO5OEAPLpKRx0i0v9UHAloVO4gphUN4dHVO4OOIiIDkIojQZ03ZTgvba2nZl9T0FFEZIBRcSSo86aOAODR1RquEpH+peJIUGXDsiktyFJxiEi/U3EkKDPjvCnDeWbDLvY1tQYdR0QGEBVHAjtv6nBa250nteihiPQjFUcCm1E8lILsdB5ZpbOrRKT/hKo4zOx8M1tnZpVmdnM3r883sxfNrM3MLgsiY5hEUoxzpwzjybU1NLe1Bx1HRAaI0BSHmUWA24ELgCnAVWY25bDdtgIfBu7q33Thdd7UETS2tPOPDbuDjiIiA0RoigOYA1S6+0Z3bwHuARZ03cHdN7v7y4DWFI85bVw+WWkRXUUuIv0mTMUxGtjW5XlVbFvczOxaM6sws4ra2uSeOE6PRjhz0jAeXbWTNt2jQ0T6QZiKo9e4+x3uXu7u5YWFhUHH6XMXTRvJ7sYWnt1YF3QUERkAwlQc1UBxl+dFsW3yNs6YOIystAh/fXl70FFEZAAIU3EsB8rMrNTM0oCFwKKAMyWEjNQI504ZzsOv7qSlTcNVItK3QlMc7t4G3AA8AqwB7nP3VWb2NTO7GMDMZptZFXA58FMzWxVc4nB5z/RR7D3YytOVu4KOIiJJLhp0gK7c/SHgocO2fbnL4+V0DmHJYeaVFTI4I8pfVm7nzEnDgo4jIkksNEcccnzSoimcf+IIHl39Ok2tuhhQRPqOiiOJXDRtFA3NbVq7SkT6lIojiZw2Lp+8rDSdXSUifUrFkUSikRQuOHEET6yp4UBLW9BxRCRJqTiSzEXTRnGwtZ0n1tQEHUVEkpSKI8nMKc1jWE46f1mp4SoR6RsqjiQTSTHePW0kT66rpf5AS9BxRCQJqTiS0GWzimhp7+DPK3TUISK9T8WRhKaOGsKUkYO5/4Vtb7+ziEicVBxJ6vLyIl6t3seaHfuCjiIiSUbFkaQWzBhNasS4v6Iq6CgikmRUHEkqLyuNcyYP508rqrViroj0KhVHEruivJi6xhaWrNU1HSLSe1QcSWxeWQHDctJ5QJPkItKLVBxJLBpJ4X0nF/G3dbXU7G8KOo6IJAkVR5K7vLyI9g7nTy/pLrwi0jtUHEluXGE2J4/J5b6KKtw96DgikgRUHAPAwtljqKxp4PlNdUFHEZEkoOIYAN4zfRSDM6L89tktQUcRkSSg4hgABqVFuLy8mMWv7tQkuYgct2jQAaR/fOCUMfx82SbufX4bnzy7LOg4kmTa2jtobGnnQEsbjc3tNLW2097htLvj7rR3QIc7He6kRVJIjaSQFu18mx7t+tzITIsSSbGgPyV5CyqOAWJsYTbzygq46/mt/OsZ40iN6GBT3lpTaztVew6yre4AVfUH2d3QzO6GFnY3NrNrfwu7GpvZ09hCY0t7r69OkJkWISs9Sk56lKz0KFnpEbLTU8lOj5CbmUZuZipDY2/zstLe9HhQagQzFU9fUnEMIB86tYRrflPB4ld38p7po4KOIyHg7tTsb2b1js4FMStrGqiqO8jWugPs3HfksObQzFTys9PJz0pj8ojBDM1KJTs9lcy0yBv/2WemRRiUGiEaMcyMiBkpZqSkgGG0dXTQ2t5BS1sHLe1OS1uX520dNLa00djcRkNzGw3N7TQ0tdLY3E51/UEamlvZe6CVfU1HvzVyWjSFobFiGZqZ1lksWankZXXmHpqVRn5W5/b8rDRyM9NIi+oXqXioOAaQsyYNo7QgizuXbeKiaSP1W9kAtL3+IMs31/Fy1V7WxMpiz4HWN14fOSSD4rxM5pYVMCYvkzF5mRTnZVI0dBB5WWmhOVJta++g/mAr9Qda2HOglbrGljce72lsYU+Xx2t27ou93nrU95eTET2iVP75OJ28w4onK21gH9WEqjjM7HzgB0AEuNPdbzns9XTgN8AsYDdwpbtv7u+ciSolxfjI6SV8+c+reHHrHmadkBd0JOlDHR1OZW0DyzfXsXxTHcs376G6/iAAGakpTByew7umjmDyyMFMHjmYiSNyGDIoNeDUPRONpFCQnU5BdnqP/86hsqlrbHnjz+7GFvZ0eVzX2Ex1fROvVO+lrrGF1vbur31Ki6Z0lkhmGvnZaW8Ml+VnpZGXnUZe7EgnP/vQMFpaUs3bhKY4zCwC3A6cC1QBy81skbuv7rLb1cAedx9vZguBbwJX9n/axHXZrCK+8+hr/HzZJhVHEtrT2MLS9bX8fV0tS9fXsquh8/bBBdnpzCkdyjXzSpldksekETlEQ3L00F/iLRt3p6G57U0FsztWMl0f725sYcvuA9Q1ttDQ3P0QWopBbmYaQ2PzMDkZqWSlR8lOj5KdHunyOEp2RvRN8zuHtmelR0MzpBaa4gDmAJXuvhHAzO4BFgBdi2MB8NXY4weA28zMXJdE91hmWpSr5ozhjqUb2LyrkZKCrKAjyXHasfcgD7+yk4de2cGLW/fQ4ZCbmcr8skLmlhUwpySPE/IzB/TQyrEwM3IyUsnJSOWE/J79nDS1tlN/oJXdjc1vOrI5/E/N/iYaamNzOM2tNLX27OSCSIqREU0hIzVC+qG3bzx+8/aMaITTxuezYMbo4/kydCtMxTEa6LqMaxVwytH2cfc2M9sL5AO7uu5kZtcC1wKMGTOmr/ImrI+eXsIvnt7ET5du4L/fNy3oOHIMduw9yEOxsnhhyx4AJo3I4YazyjhzYiHTinKTamgkUWSkRhgxJMKIIRlx/b1DpzM3NHc5MaCp60kCnY8PtrbT3NpBU1s7Ta0dNLd10NTaefpzc2sHexpbaHrj9XbystNY0AefZ5iKo9e4+x3AHQDl5eU6GjnMsMEZXFFexL3Lt/GpsyfE/U0uwTjQ0saDL+/gvoptLN/cWRaTRw7mc+dN4MKTRjK2MDvghHKsopEUhgxKSZw5pqADdFENFHd5XhTb1t0+VWYWBYbQOUkucbpu/jjufn4bP3tqI//voilBx5GjcHdWVu3l3uXb+MvK7TQ0tzG2MEtlIYEKU3EsB8rMrJTOglgIvP+wfRYBHwKeAS4Dlmh+49gU52WyYPoo7npuKx8/Y1xcZ6dI32tobuOBim3cs3wba3fuZ1BqhHdPG8mVs4spP2Go5iskUKEpjticxQ3AI3SejvsLd19lZl8DKtx9EfBz4LdmVgnU0VkucoyuP2s8f1pRzU+e3MC/66gjFKr2HODX/9jMPc9vY39zG9OKhvCNS07iPdNHkpORGMMYkvxCUxwA7v4Q8NBh277c5XETcHl/50pW4wqzed/JRfz22S1cM2+s5joC9NLWPdy5bBOLX90JwIUnjeTquaXMKM4NNphIN0JVHNL/PnV2GX96qZrb/1bJ1997YtBxBhR358nXarltSSUvbNlDTkaUa+aW8sHTShidOyjoeCJHpeIY4IrzMrlydjH3LN/KNfNKe3y+uhw7d+eJNTX8cMl6Xq7ay+jcQXz1PVO4rLyY7HT9SEr46btU+NTZZfzxpWq+uXgtP/rArKDjJK2ODufR1a/zP0vWs2r7PorzBvHNS0/ikplFobkiWKQnVBzCsMEZXDd/HN97/DWWb65jdomWIulN7s6StTV8+5F1rN25n5L8TG69fDoLZowKzaKBIvHQd60A8LH5pYwYnMF//nU1HR06w7m3LN9cx+U/eYarf11Bc1sH379yBo9/5p1cNqtIpSEJS0ccAnSuYfX5d03ks/ev5P4XtnHlbC3VcjzW7NjHrY+s44m1NQzLSecbl5zE5eUqC0kOKg55wyUzR3PP8q3898NrOXfKCPKy0oKOlHCq6w9y6yPr+NOKanLSo9x0/iQ+fFoJg9IiQUcT6TX69UfekJJi/Od7T6KhqY1vPrw26DgJpbG5je88uo6zbn2Sh17ZwXXzx/HUF87i42eMU2lI0tERh7zJxBE5XD23lJ8u3ciCmaM4bVxB0JFCraPD+cNL1Xxr8Vpq9jdz8fRR3HTBJF2HIUlNRxxyhE+fM4HSgiy+8MDLR70xjXROfL/3R0/zuftXMjJ3EP/78VP54VUzVRqS9FQccoRBaRFuvXwa2+sP8l8Prgk6TuhsqzvA9Xe9yOU/eYaafc1878rp/PHjp+mOijJgaKhKujXrhDw+Nm8sP126kbMmDePcKcODjhS4huY2fvxkJT97ahMp1nnh5HXvHEtmmn6MZGDRd7wc1WfOm8DTG3bxuftX8uCNcykamhl0pEB0dDgPvFjFtx9ZR+3+Zt47YxRfOH8SozQkJQOUhqrkqNKjEW676mTaO5wb736J1vae3Rc5mTy3cTcX376MLzzwMqNzB/GHT5zG9xfOVGnIgKbikLdUUpDFLZeexItb6/nKolUMlPtmbd19gI//7gWuvONZdje08IOFM/jjJ07j5DFDg44mEjgNVcnbumjaKFZt38ePn9zApBE5fPDUkqAj9Zn9Ta3c9rdKfrlsM5EU4zPnTuBj88bqWgyRLlQc0iOfP28i61/fz3/8ZTWjcwdx9uTkmixv73Duq9jGdx5dx66GFi49uYjPv2uibm4l0g0NVUmPpKQY3184k6mjBvOJ37/Icxt3Bx2p1zxduYt3//ApvviHVygtyGLRDafznSumqzREjkLFIT2WnR7lVx+ZQ9HQQVz96wpe2ron6EjHZfX2fXz0V8v5wJ3P0dDcxo8+cDL3XXcq04pyg44mEmoqDolLXlYav736FPKy0viXO5/j+U11QUeK2+Zdjdx490tc+MOnqNhcx80XTOLxz7yTC08aiZkFHU8k9CzZz5IpLy/3ioqKoGMknZ17m3j/nc+yo76J294/MyHmPF7f18QPn1jPvcu3kRpJ4SOnl3Dd/HEMyUwNOppI6JjZC+5e3t1rmhyXYzJiSAb3XnsqH/nV81zzmwq+dOFkrp5bGsrf2Gv2NfGzpzby22e30NbuXDVnDJ88azzDBmsOQ+RYqDjkmBXmpHP/dafxmftW8J8PrmHFtnq+8b6TGJwRjt/gt+4+wE+XbuD+iiraOjpYMGM0nz6njBPys4KOJpLQQlEcZpYH3AuUAJuBK9z9iJlXM1sMvANY5u4X9WdG6d6gtAi3v/9kfvz3DXz3sdd4uWov375sGqeMzQ8kj7tTsWUPv3lmCw+9soOIGZeVF3Hd/LEqDJFeEoo5DjP7FlDn7reY2c3AUHe/qZv9zgYyget6Whya4+g/FZvr+PS9K6jac5CFs4v5/Lsmkp+d3i8fu6G5jT+9VM3vnt3C2p37GZwR5crZxVwzbyzDNSQlEre3muMIS3GsA85w9x1mNhJ40t0nHmXfM4DPqTjC6UBLG99/fD0/X7aJQakRrplXykdOL2XIoN4fvmpua+fJdbUsWrmdJ9a8TlNrB1NGDuaDp57AghmjdbW3yHFIhOKod/fc2GMD9hx63s2+Z/A2xWFm1wLXAowZM2bWli1bejmxvJ3Kmv3c+shrLF61k0GpES6dNZrLZxUzrWjIcU2gb68/yNOVu3i6chdPrK1hf1MbeVlpXHjSCC6ZWcTJY3JDOUEvkmhCURxm9jgwopuXvgT8umtRmNked+92NTkdcSSWVdv38sunN7NoxXZa2jsoGjqIeWWFzCkdypSRQzghP5OM1COPDDo6nO17D7K+poENNQ289vp+KjbvYeOuRgDys9J458RCLp4+itPHF5Aa0SVJIr0pFMXxVjRUlfz2Hmjl0dU7eWTVTp7bVMf+ps5b0ppB7qBUhgxKJZJidDjsb2qj/kALbR3//N7Mz0pjWtEQTh9fwOnjC5g4PIeUFB1ZiPSVRLiOYxHwIeCW2Ns/BxtHetuQzFQuLy/m8vJi2juc117fz/qaBjbWNrC7oYW9B1tpdyfFjOz0KLmZqRQNHUTZsBzGD8smLyst6E9BRGLCUhy3APeZ2dXAFuAKADMrB/7V3a+JPX8KmARkm1kVcLW7PxJQZjlGkRRj8sjBTB45OOgoInIMQlEc7r4bOLub7RXANV2ez+vPXCIiciTNKIqISFxUHCIiEhcVh4iIxEXFISIicVFxiIhIXFQcIiISFxWHiIjEJRRLjvQlM6ul86LCY1EA7OrFOH1BGY9f2PNB+DOGPR8oY7xOcPfC7l5I+uI4HmZWcbS1WsJCGY9f2PNB+DOGPR8oY2/SUJWIiMRFxSEiInFRcby1O4IO0APKePzCng/CnzHs+UAZe43mOEREJC464hARkbioOEREJC4qjqMws/PNbJ2ZVZrZzUHnOZyZFZvZ38xstZmtMrNPBZ2pO2YWMbOXzOyvQWfpjpnlmtkDZrbWzNaY2alBZ+rKzP5v7N/3VTO728wyQpDpF2ZWY2avdtmWZ2aPmdn62NuhIcz47di/88tm9kczyw1Tvi6vfdbM3MwKgsjWEyqObphZBLgduACYAlxlZlOCTXWENuCz7j4FeAdwfQgzAnwKWBN0iLfwA2Cxu08CphOirGY2GrgRKHf3E4EIsDDYVAD8Cjj/sG03A0+4exnwROx5kH7FkRkfA05092nAa8AX+ztUF7/iyHyYWTFwHrC1vwPFQ8XRvTlApbtvdPcW4B5gQcCZ3sTdd7j7i7HH++n8D290sKnezMyKgHcDdwadpTtmNgSYD/wcwN1b3L0+0FBHigKDzCwKZALbA86Duy8F6g7bvAD4dezxr4H39memw3WX0d0fdfe22NNngaJ+D/bPLN19DQG+B3wBCPVZSyqO7o0GtnV5XkXI/lPuysxKgJnAcwFHOdz36fwh6Ag4x9GUArXAL2PDaXeaWVbQoQ5x92rgVjp/+9wB7HX3R4NNdVTD3X1H7PFOYHiQYXrgo8DDQYfoyswWANXuvjLoLG9HxZHgzCwb+F/g0+6+L+g8h5jZRUCNu78QdJa3EAVOBn7s7jOBRoIfYnlDbJ5gAZ0FNwrIMrN/CTbV2/POc/xD+xuzmX2JzqHe3wed5RAzywT+Dfhy0Fl6QsXRvWqguMvzoti2UDGzVDpL4/fu/oeg8xzmdOBiM9tM51DfWWb2u2AjHaEKqHL3Q0dqD9BZJGFxDrDJ3WvdvRX4A3BawJmO5nUzGwkQe1sTcJ5umdmHgYuAD3i4LmIbR+cvCCtjPzNFwItmNiLQVEeh4ujecqDMzErNLI3OCclFAWd6EzMzOsfm17j7d4POczh3/6K7F7l7CZ1fvyXuHqrflt19J7DNzCbGNp0NrA4w0uG2Au8ws8zYv/fZhGjy/jCLgA/FHn8I+HOAWbplZufTOXR6sbsfCDpPV+7+irsPc/eS2M9MFXBy7Hs0dFQc3YhNoN0APELnD+p97r4q2FRHOB34P3T+Jr8i9ufCoEMloE8Cvzezl4EZwDeCjfNPsSOhB4AXgVfo/HkNfEkKM7sbeAaYaGZVZnY1cAtwrpmtp/NI6ZYQZrwNyAEei/28/CRk+RKGlhwREZG46IhDRETiouIQEZG4qDhERCQuKg4REYmLikNEROKi4hARkbioOEREJC4qDpF+YGZFZnblYdtKzOyh2H1fXjOzIJf5FukxFYdI/zibLutgmVkKneuM/cTdJwInAeVmdm1A+UR6TFeOi/QxM5tL59pN9cB+4H3AROAad7+0y34jgb+7+4Qgcor0lI44RPqYuy+jc+HMBe4+w903ApOBlYfttwMYHFtYUyS0VBwi/WMisLbL83Ygu+sOsRVwM+m8V4RIaKk4RPqYmRXQefe+roXwJHBhrCwOORd40d3DesdEEUDFIdIfSjjsXuGx24O+BHwNwMyGA9+l8y5wIqGm4hDpe2uBAjN71cxOAzCzm4Fy4N/N7Czgx8AJwI9i95AXCS2dVSUiInHREYeIiMRFxSEiInFRcYiISFxUHCIiEhcVh4iIxEXFISIicVFxiIhIXP4/anJn+FVX4eMAAAAASUVORK5CYII=\n", "text/plain": [ "
" ] }, "metadata": { "needs_background": "light" }, "output_type": "display_data" } ], "source": [ "Omega = 1.0\n", "omega_cutoff = 5.0\n", "alpha = 0.3\n", "\n", "system = oqupy.System(0.5 * Omega * sigma_x)\n", "correlations = oqupy.PowerLawSD(alpha=alpha,\n", " zeta=1,\n", " cutoff=omega_cutoff,\n", " cutoff_type='exponential')\n", "bath = oqupy.Bath(0.5 * sigma_z, correlations)\n", "tempo_parameters = oqupy.TempoParameters(dt=0.1, tcut=3.0, epsrel=10**(-4))\n", "\n", "dynamics = oqupy.tempo_compute(system=system,\n", " bath=bath,\n", " initial_state=up_density_matrix,\n", " start_time=0.0,\n", " end_time=15.0,\n", " parameters=tempo_parameters)\n", "t, s_z = dynamics.expectations(0.5*sigma_z, real=True)\n", "\n", "plt.plot(t, s_z, label=r'$\\alpha=0.3$')\n", "plt.xlabel(r'$t\\,\\Omega$')\n", "plt.ylabel(r'$\\langle\\sigma_z\\rangle$')\n", "plt.legend()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### 1. The model and its parameters \n", "We consider a system Hamiltonian\n", "$$ H_{S} = \\frac{\\Omega}{2} \\hat{\\sigma}_x \\mathrm{,}$$\n", "a bath Hamiltonian\n", "$$ H_{B} = \\sum_k \\omega_k \\hat{b}^\\dagger_k \\hat{b}_k \\mathrm{,}$$\n", "and an interaction Hamiltonian\n", "$$ H_{I} = \\frac{1}{2} \\hat{\\sigma}_z \\sum_k \\left( g_k \\hat{b}^\\dagger_k + g^*_k \\hat{b}_k \\right) \\mathrm{,}$$\n", "where $\\hat{\\sigma}_i$ are the Pauli operators, and the $g_k$ and $\\omega_k$ are such that the spectral density $J(\\omega)$ is\n", "$$ J(\\omega) = \\sum_k |g_k|^2 \\delta(\\omega - \\omega_k) = 2 \\, \\alpha \\, \\omega \\, \\exp\\left(-\\frac{\\omega}{\\omega_\\mathrm{cutoff}}\\right) \\mathrm{.} $$\n", "Also, let's assume the initial density matrix of the spin is the up state\n", "$$ \\rho(0) = \\begin{pmatrix} 1 & 0 \\\\ 0 & 0 \\end{pmatrix} $$\n", "and the bath is initially at zero temperature." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "For the numerical simulation it is advisable to choose a characteristic frequency and express all other physical parameters in terms of this frequency. Here, we choose $\\Omega$ for this and write:\n", " \n", "* $\\Omega = 1.0 \\Omega$\n", "* $\\omega_c = 5.0 \\Omega$\n", "* $\\alpha = 0.3$" ] }, { "cell_type": "code", "execution_count": 5, "metadata": {}, "outputs": [], "source": [ "Omega = 1.0\n", "omega_cutoff = 5.0\n", "alpha = 0.3" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### 2. Create system, correlations and bath objects" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "#### System\n", "$$ H_{S} = \\frac{\\Omega}{2} \\hat{\\sigma}_x \\mathrm{,}$$" ] }, { "cell_type": "code", "execution_count": 6, "metadata": {}, "outputs": [], "source": [ "system = oqupy.System(0.5 * Omega * sigma_x)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "#### Correlations\n", "$$ J(\\omega) = 2 \\, \\alpha \\, \\omega \\, \\exp\\left(-\\frac{\\omega}{\\omega_\\mathrm{cutoff}}\\right) $$" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Because the spectral density is of the standard power-law form,\n", "$$ J(\\omega) = 2 \\alpha \\frac{\\omega^\\zeta}{\\omega_c^{\\zeta-1}} X(\\omega,\\omega_c) $$\n", "with $\\zeta=1$ and $X$ of the type ``'exponential'`` we define the spectral density with:" ] }, { "cell_type": "code", "execution_count": 7, "metadata": {}, "outputs": [], "source": [ "correlations = oqupy.PowerLawSD(alpha=alpha,\n", " zeta=1,\n", " cutoff=omega_cutoff,\n", " cutoff_type='exponential')" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "#### Bath\n", "The bath couples with the operator $\\frac{1}{2}\\hat{\\sigma}_z$ to the system." ] }, { "cell_type": "code", "execution_count": 8, "metadata": {}, "outputs": [], "source": [ "bath = oqupy.Bath(0.5 * sigma_z, correlations)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### 3. TEMPO computation\n", "Now, that we have the system and the bath objects ready we can compute the dynamics of the spin starting in the up state, from time $t=0$ to $t=5\\,\\Omega^{-1}$" ] }, { "cell_type": "code", "execution_count": 9, "metadata": {}, "outputs": [ { "name": "stderr", "output_type": "stream", "text": [ "../oqupy/tempo.py:881: UserWarning: Estimating parameters for TEMPO computation. No guarantee that resulting TEMPO computation converges towards the correct dynamics! Please refer to the TEMPO documentation and check convergence by varying the parameters for TEMPO manually.\n", " warnings.warn(GUESS_WARNING_MSG, UserWarning)\n", "WARNING: Estimating parameters for TEMPO computation. No guarantee that resulting TEMPO computation converges towards the correct dynamics! Please refer to the TEMPO documentation and check convergence by varying the parameters for TEMPO manually.\n" ] }, { "name": "stdout", "output_type": "stream", "text": [ "--> TEMPO computation:\n", "100.0% 80 of 80 [########################################] 00:00:01\n", "Elapsed time: 2.0s\n" ] } ], "source": [ "dynamics_1 = oqupy.tempo_compute(system=system,\n", " bath=bath,\n", " initial_state=up_density_matrix,\n", " start_time=0.0,\n", " end_time=5.0,\n", " tolerance=0.01)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "and plot the result:" ] }, { "cell_type": "code", "execution_count": 10, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "" ] }, "execution_count": 10, "metadata": {}, "output_type": "execute_result" }, { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAY4AAAEICAYAAABI7RO5AAAAOXRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjMuNCwgaHR0cHM6Ly9tYXRwbG90bGliLm9yZy8QVMy6AAAACXBIWXMAAAsTAAALEwEAmpwYAAAnjklEQVR4nO3dd3yV5f3/8dcnexPIYCRA2EOmnDBEcYEiWpW6AAVUEK1QtVOt7beOn98O+7XWDYKKiiK2VWlLRawDFyMgyhKIDAkzhBVG9vX7I9FGDCOQk/vknPfz8cgjuc+5cvI+iLxz3/d1X7c55xARETlRYV4HEBGRhkXFISIitaLiEBGRWlFxiIhIrag4RESkViK8DuBvqampLisry+sYIiINypIlS3Y559Jqei7oiyMrK4ucnByvY4iINChmtuloz+lQlYiI1IqKQ0REakXFISIitRL05zhEREpLS8nLy6OoqMjrKAEnJiaGzMxMIiMjT/h7VBwiEvTy8vJITEwkKysLM/M6TsBwzlFQUEBeXh5t2rQ54e8LqENVZjbUzNaYWa6Z3VXD89ebWb6ZLav6GO9FThFpWIqKikhJSVFpHMHMSElJqfWeWMDscZhZOPAEMATIAxab2Wzn3Kojhr7qnJtU7wFFpEFTadTsZP5cAqY4gL5ArnNuPYCZzQQuA44sjnqxqeAgf1uSR+P4KJpUfTRvFEvb1HjCwvQXUERCVyAVRwawudp2HtCvhnFXmNkgYC3wE+fc5iMHmNkEYAJAq1atTirM+vyDPPZeLkferiQ5LhJf68ZkZzXhrA5pdG2RdFKvLyLSUAVScZyIfwCvOOeKzexmYDpw3pGDnHNTgCkAPp/vpO5UdW7ndHIfHMa+w6XsPljM7oOlbCw4SM7G3SzeuId3Vu/kd//+ktNaJHG1ryWX9WpBclzUqbw3EZEGIZCKYwvQstp2ZtVj33LOFVTbnAr80Z+BwsPs28NUAH3bNOFqX2XEnYVFvLViO7NyNvPb2St58F+ruaxXC247vwMtm8T5M5aIhLi33nqL22+/nfLycsaPH89dd31vLhFFRUUMGjSI4uJiysrKuPLKK7nvvvvq5OcHUnEsBjqYWRsqC2MEMKr6ADNr7pzbVrV5KbC6fiP+V3piDGMGZDFmQBYrt+5j5qLNvJqzmTeWbeFqX0smndee5o1ivYonIkGqvLyciRMnMm/ePDIzM8nOzubSSy+la9eu3xkXHR3Nu+++S0JCAqWlpZx55plcdNFF9O/f/5QzBMx0XOdcGTAJmEtlIcxyzq00s/vN7NKqYbeZ2Uoz+xy4Dbjem7TfdVqLRjxweTc++MU5XJPdklk5mzn7ofd5+O01FJWWex1PRALEypUrGTx4MB07duSBBx7gxz/+MYsXL67VayxatIj27dvTtm1boqKiGDFiBG+++eb3xpkZCQkJQOUFkKWlpXU2syyQ9jhwzs0B5hzx2P9U+/pu4O76znWimjeK5f9d3p2bB7XjT2+v4dF3c/nnF9t4cHh3BrRL8TqeiAD3/WMlq7bur9PX7Noiid/+4LRjjikqKuKqq67itddeo23btnTu3Jk+ffqQnZ397ZizzjqLwsLC733vn/70JwYPHgzAli1baNnyv0f1MzMzWbhwYY0/s7y8nD59+pCbm8vEiRPp16+m+Ua1F1DFESxaNonjLyN6c2WfTO55fQUjn1nA1b5Mfn1JV5JiTvyyfhEJHu+88w69e/fmtNMqC6akpISf/exn3xnz4Ycf1unPDA8PZ9myZezdu5fhw4ezYsUKunXrdsqvq+Lwo7M6pDH3jkE8+u46psxfz6frC3h0RG96t2rsdTSRkHW8PQN/WbZsGb179wZg69atJCQkMHDgwO+MOZE9joyMDDZv/u9VCHl5eWRkZBzzZycnJ3Puuefy1ltvqTgagtiocO4c2pnBXZpy2yufcdXTn/KzCzpx86C2upBQJIRERUWxZUvlRNG7776bkpKS7405kT2O7Oxs1q1bx4YNG8jIyGDmzJm8/PLL3xuXn59PZGQkycnJHD58mHnz5nHnnXee+hshgE6OB7s+rRsz5/azuPC0ZvzhrS8Z+9wi9hz8/l8cEQlOo0aNYv78+XTq1ImePXsyYMAA7rjjjlq/TkREBI8//jgXXnghXbp04eqrr/728BfAsGHD2Lp1K9u2bePcc8+lR48eZGdnM2TIEC655JI6eS/mjrw0Osj4fD4XSLeOdc7xyqLN3Dt7Jc0axTB1rI+OTRO9jiUS1FavXk2XLl28jhGwavrzMbMlzjlfTeO1x1HPzIxR/Vox8+b+HC4tZ/gTHzNv1Q6vY4mInDAVh0dOb9WYf0w6k3bpCUx4MYcn3ssl2Pf+RCQ4qDg81KxRDLNuHsAPerTgoblr+PUbKyivUHmI+IN+MavZyfy5aFaVx2Iiw3nkml60SI7l6Q++Ir+wmEdH9iYmMtzraCJBIyYmhoKCAt3M6Qjf3AEwJiamVt+n4ggAYWHGXRd1pllSNPf9cxXXTV3I1LE+rbYrUkcyMzPJy8sjPz/f6ygB55t7jteGiiOAXD+wDWmJMfzk1WVcM3kBL47vS3pi7X4TEJHvi4yMrNU9teXYdI4jwFzcoznP35DN17sPMWLyArbtO+x1JBGR71BxBKAz2qfy4ri+5BcWc9XTn7J59yGvI4mIfEvFEaB8WU2YcVM/DhSXcdXTn7I+/4DXkUREABVHQOuRmczMCf0pq6hg5DML2LDroNeRRERUHIGuc7MkXr6pP6XljpFTFrCpQOUhIt5ScTQAHZsmMmN8P4rLyhk5ZYHOeYiIp1QcDUSX5km8NL4fB0vKGTFlAVv2araViHhDxdGAnNaiES+N68f+olKum7qQ/MJiryOJSAhScTQw3TMb8fwN2WzfV8ToaQvZd6jU60giEmJUHA1Qn9ZNmDKmD+vzD3L984s4WFzmdSQRCSEqjgbqrA5pPDqyN1/k7eOmF3IoLiv3OpKIhAgVRwM2tFszHrqyB598VcBPXl2mJdlFpF4EVHGY2VAzW2NmuWZ21zHGXWFmzsxqvK1hKPnh6Zn8+uIuzFm+nd/OXqF7DoiI3wXM6rhmFg48AQwB8oDFZjbbObfqiHGJwO3AwvpPGZjGn9WW/APFTP5gPWkJMdw+uIPXkUQkiAXSHkdfINc5t945VwLMBC6rYdwDwB+AovoMF+juGtqZK07P5M/vrOWlBZu8jiMiQSyQiiMD2FxtO6/qsW+Z2elAS+fcv471QmY2wcxyzCwnVG7cYmb8/orunNc5nf95cwVvr9zudSQRCVKBVBzHZGZhwMPAz4431jk3xTnnc8750tLS/B8uQESGh/H4qN50z2jEbTM/47Ov93gdSUSCUCAVxxagZbXtzKrHvpEIdAPeN7ONQH9gtk6Qf1dcVATTrs8mPTGGcdNz2KgVdUWkjgVScSwGOphZGzOLAkYAs7950jm3zzmX6pzLcs5lAQuAS51zOd7EDVypCdFMv7EvzjnGPreIggNamkRE6k7AFIdzrgyYBMwFVgOznHMrzex+M7vU23QNT5vUeKZdX7k0ybjpORSV6gJBEakbFuzz/n0+n8vJCd2dkrdWbOdHM5Yw9LRmPDHqdMLCzOtIItIAmNkS51yNpwICZo9D/GNot2bcM6wL/16xnT/M/dLrOCISBALmAkDxn3FntmFTwSEmf7Ce1k3iGdWvldeRRKQBU3GEADPjtz/oSt6eQ/zmzRVkNI7l7I6hM01ZROqWDlWFiIjwMB4bdTodmyYyacZS1u4o9DqSiDRQKo4QkhAdwbSxPmKiwrnx+cXs0jRdETkJKo4Q0yI5lqljfOw6UMyEFzRNV0RqT8URgnq2TObhq3ux9Ou93Pm3L7QUu4jUioojRA3r3pxfXNiJN5dt5bF3c72OIyINiGZVhbBbz2nHV/kHeHjeWtqnJzCse3OvI4lIA6A9jhBmZvzv8O6c3iqZn85axvK8fV5HEpEGQMUR4mIiw5k82kdKfDQ3vZDDjv26P5aIHJuKQ0hLjGbqWB/7i0o100pEjkvFIQB0aZ7EX0b05ost+/jlXzXTSkSOTsUh3xrStSm/uLATsz/fypPvf+V1HBEJUJpVJd/xo7PbsXZ7IQ/NXUOH9AQuOK2Z15FEJMBoj0O+w8z4/RU96JnZiDteXcaX2/d7HUlEAoyKQ74nJjKcKWN8JMZEMH56jm49KyLfoeKQGjVNimHKaB/5hcX8aMZSSsoqvI4kIgFCxSFH1bNlMn+8sgeLNuzm3n+s1EwrEQF0clyO47JeGXy5vZCn3v+KLs0SGT0gy+tIIuIx7XHIcf38gk6c3zmde/+xik9yd3kdR0Q8puKQ4woPMx4Z0Yu2qfHc+vJSvi445HUkEfGQikNOSGJMJM+M8eEc3PRCDgeKy7yOJCIeCajiMLOhZrbGzHLN7K4anr/FzJab2TIz+8jMunqRM1Rlpcbz+KjerNtZyE9fXUZFhU6Wi4SigCkOMwsHngAuAroCI2sohpedc92dc72APwIP129KOatDGvdc3JW3V+3gkf+s8zqOiHggYIoD6AvkOufWO+dKgJnAZdUHOOeqX8YcD+hXXg/cODCLK/tk8uh/1jFn+Tav44hIPQuk4sgANlfbzqt67DvMbKKZfUXlHsdtNb2QmU0wsxwzy8nPz/dL2FBmZjw4vBu9WyXzs1mfs3KrbgAlEkoCqThOiHPuCedcO+BO4NdHGTPFOedzzvnS0tLqN2CIiI4IZ/J1fWgUG8mEF5ZoWRKREBJIxbEFaFltO7PqsaOZCVzuz0BybOlJMUwe3Yf8A1qWRCSUBFJxLAY6mFkbM4sCRgCzqw8wsw7VNi8GdHbWYz1bJvPHKyqXJbnvHyu9jiMi9SBglhxxzpWZ2SRgLhAOPOucW2lm9wM5zrnZwCQzGwyUAnuAsd4llm9c3juD1dv2M3n+ero0T+K6/q29jiQifhQwxQHgnJsDzDnisf+p9vXt9R5KTsgvh3Zm7Y5C7p29kvbpCfRvm+J1JBHxk0A6VCUNWHiY8ZeRvWmVEsetM5ayebeWJREJVioOqTNJMZFMHeOjtLyCm17I4aCWJREJSioOqVNt0xJ4fNTprN1RyM9f+1zLkogEIRWH1LmzO6bxq2Fd+PeK7Tz6ria+iQSbgDo5LsFj3JltWL2tkEfeWUfnZokM7dbc60giUke0xyF+UX1Zkp+8+jmrtu4//jeJSIOg4hC/iYn877IkN72Qo2VJRIKEikP8Kj0philj+rBLy5KIBA0Vh/hdj8xk/nhl5bIkv529Euc000qkIdPJcakXl/XKYM32Qp58/ys6N0tk7BlZXkcSkZOkPQ6pNz+/oBODuzTl/n+u4uPcXV7HEZGTpOKQehMWZjwyohft0uK5dcZSNuw66HUkETkJKg6pVwnREUwdk02Ywfjpi9l3uNTrSCJSSyoOqXetUuJ46ro+bCo4xI9f+Yyycs20EmlIVBziif5tU3jg8m7MX5vP/8750us4IlILmlUlnhnZtxVrdxTy7Mcb6Ng0gRF9W3kdSUROgPY4xFP3DOvCoI5p/ObNFSxcX+B1HBE5ASoO8VREeBiPjexNqyZx3PLSEjYVaKaVSKBTcYjnGsVGMm1sNg4YNz2H/UWaaSUSyFQcEhCyUuN56to+bNx1kEkva6aVSCBTcUjAGNAuhQeHV860+n//Wu11HBE5Cs2qkoByTXYr1u04wNSPNtAuLZ7RA7K8jiQiR1BxSMC5e1gXNuw6yL3/WEWrlHjO7pjmdSQRqeaEDlWZWWMzG2Nmr5vZSjP7p5ndZGbpdRnGzIaa2RozyzWzu2p4/qdmtsrMvjCz/5hZ67r8+RIYwsOMv4zsTYf0BCbNWMraHYVeRxKRao5bHGb2d+DvQDpwp3PuNOBWIAZ4yczer4sgZhYOPAFcBHQFRppZ1yOGfQb4nHM9gL8Cf6yLny2BJyE6gmevzyYmKpwbn1/MLt09UCRgnMgex43OuXOdc39yzq0FcM597Zx7zDl3AXB5HWXpC+Q659Y750qAmcBl1Qc4595zzh2q2lwAZNbRz5YA1CI5lmfG+MgvLGbCCzkUlZZ7HUlEOIHicM7tPfIxM3vRzGLqOEsGsLnadl7VY0czDvh3TU+Y2QQzyzGznPz8/DqMKPWtV8tk/nxNL5Z+vZefv/Y5FRW6e6CI1052Oq4BT5tZLPDTOsxzYj/c7DrABzxU0/POuSnOOZ9zzpeWphOrDd2w7s25c2hn/vnFNh6et9brOCIh72RnVW0ApgFPAXW1wNAWoGW17cyqx77DzAYD9wBnO+d04DtE3HJ2WzYVHOTx93JpnRLHVb6Wx/8mEfGLk93jmOKc2wjcCwytoyyLgQ5m1sbMooARwOzqA8ysNzAZuNQ5t7OOfq40AGbGA5d348z2qdz99+V8olvPinjmhIvDzO795mvn3OaqzxurZlmdMudcGTAJmAusBmY551aa2f1mdmnVsIeABOA1M1tmZrOP8nIShCLDw3jyutNpkxrPLS8tYZ2m6Yp4wpw7sZONZlZB5T/cTYClwEzn3B4/ZqsTPp/P5eTkeB1D6lDenkMMf/ITosLDeH3iGaQn1vU8DRExsyXOOV9Nz9XmUJUDiqjcI2gJfGJmPesgn0itZDaO49mx2ew+WMK453M4VFLmdSSRkFKb4vjSOfdb59xfnXO/ovIaiz/7KZfIMXXPbMTjo3qzcus+bnvlM8o1TVek3tSmOHaZWZ9vNqouBtRcV/HM+V2act9l3Xhn9U7unb2SEz3sKiKnpjbTcW8DZprZEmA50JPKabkinhndvzV5ew4x+YP1tEiO5UfntPM6kkjQO25xmFmYc67COfe5mfUCBgPdgHeBl7953s85RY7qzgs7s21vEX9460uaN4rh8t7HWnBARE7ViRyqmmdmr5rZSCDaOfcv4DFgD5XXVCz1Z0CR4wkLMx66qgcD2qbwi79+zse6xkPEr05krarzgfuALOBfZrYA+A/QHfizc66XPwOKnIjoiHCeHt2HtqkJ3PziElZt3e91JJGgdcLXcXz7DWaxzrnDfspT53QdR2jZtu8wP3zyE8orHH/70Rm0bBLndSSRBqmuruMAoCGVhoSe5o1imX5jX4pKyxn77CJ2HyzxOpJI0DnZtapEAlbHpolMuz6bLXsPc8Pzi3WBoEgdU3FIUMrOasJjI3uzPG8vt85YSmm5Jv6J1BUVhwStC05rxoPDu/P+mnzu/OsXugmUSB052ftxiDQII/u2YldhMf83by2N46P49cVdMDOvY4k0aCoOCXqTzmtPwcESpn20gSbxUUw8t73XkUQaNBWHBD0z438u6creQyU8NHcNTeKjGNm3ldexRBosFYeEhMqry3uy93Ap97y+nKSYSC7u0dzrWCINkk6OS8iIDA/jqWv70Kd1Y+549TPeX6O7D4ucDBWHhJTYqHCmXZ9Nx6aJ3PLSEhZt2O11JJEGR8UhIScpJpLpN/alRXIs455fzIot+7yOJNKgqDgkJKUmRPPSuH4kxUYy5tlFrNtR6HUkkQZDxSEhq0VyLDPG9yM8zLh26kI27jrodSSRBkHFISEtKzWel8f3o6zCce3UheTtOeR1JJGAF1DFYWZDzWyNmeWa2V01PD/IzJaaWZmZXelFRgk+HZom8sKNfSksKuXaqQvZsb/I60giAS1gisPMwoEngIuArsBIM+t6xLCvgeuBl+s3nQS7bhmNmH5jX3YVFjPqmQXkFxZ7HUkkYAVMcQB9gVzn3HrnXAkwE7is+gDn3Ebn3BeAljqVOte7VWOeu6EvW/cWMeqZBew6oPIQqUkgFUcGsLnadl7VY7VmZhPMLMfMcvLz8+sknISGvm2a8Oz12Wzec4jrpi7UjaBEahBIxVFnnHNTnHM+55wvLS3N6zjSwAxol8K0sdls2HWQa6cuZO8hlYdIdYFUHFuAltW2M6seE6l3A9un8swYH1/lH1B5iBwhkIpjMdDBzNqYWRQwApjtcSYJYYM6pjFldB/W7TzAqGcWskeHrUSAACoO51wZMAmYC6wGZjnnVprZ/WZ2KYCZZZtZHnAVMNnMVnqXWELBOZ3SeWaMj9z8A4zSOQ8RAMy54L6dps/nczk5OV7HkAZu/tp8bnohhzap8cwY34+UhGivI4n4lZktcc75anouYPY4RALZoI5pTB3rY8Oug4yYsoCdhbpIUEKXikPkBJ3VIY3nbshmy97DjJi8gG37DnsdScQTKg6RWjijXSov3NiXnYXFXD35Uzbv1tpWEnpUHCK15Mtqwozx/dh/uIyrJ3/K+vwDXkcSqVcqDpGT0LNlMq/c1J+Ssgqunvwpq7bu9zqSSL1RcYicpK4tkph1ywAiw8MYMeVTlmza43UkkXqh4hA5Be3SEnjtlgE0iY/iuqkL+WjdLq8jifidikPkFGU2jmPWLQNonRLHjc8vZs7ybV5HEvErFYdIHUhPjOHVCQPokdmIiS8v5aUFm7yOJOI3Kg6ROtIoLpIXx/XjvE7p/PqNFTzyzlqCfWUGCU0qDpE6FBsVztOj+3DF6Zk88s46fvPmCsorVB4SXCK8DiASbCLDw/jTVT1ITYhi8vz17NhfzKMjehMbFe51NJE6oT0OET8wM+4e1oV7f9CVd1bvYNTUBVpZV4KGikPEj64f2Ianru3Dqq37ueKpT9hUcNDrSCKnTMUh4mdDuzXj5Zv6sedQCcOf/IQlm3Z7HUnklKg4ROpBn9ZNeP3WgTSKjWTkMwt5c5nuiiwNl4pDpJ60SY3n7z86g14tk7l95jIe/c86TdeVBknFIVKPGsdH8eK4vvzw9AwenreW22cuo6i03OtYIrWi6bgi9Sw6Ipz/u6on7dIS+NPba9hYcJApo300axTjdTSRE6I9DhEPmBkTz23PlNE+vtp5gB88/hGffa3VdaVhUHGIeGhI16b8/daBxESGcc2UBcxavNnrSCLHpeIQ8VinZom8OfFMsrMa88u/fcGvXl9OcZnOe0jgUnGIBIAm8VFMv6EvN5/dlpcXfs2IKQvYvq/I61giNQqo4jCzoWa2xsxyzeyuGp6PNrNXq55faGZZHsQU8YuI8DDuvqgLT157Omu3F3LJYx/yca5uDCWBJ2CKw8zCgSeAi4CuwEgz63rEsHHAHudce+DPwB/qN6WI/w3r3pw3Jg4kOS6K66Yt5JF31mqFXQkoAVMcQF8g1zm33jlXAswELjtizGXA9Kqv/wqcb2ZWjxlF6kWHponMnjSQ4b0yeOSddYx5diH5hcVexxIBAqs4MoDqU0ryqh6rcYxzrgzYB6Qc+UJmNsHMcswsJz8/309xRfwrLiqC/7u6J3+4ojs5G/dw0V8+5MN1+vss3guk4qgzzrkpzjmfc86XlpbmdRyRk2ZmXJPdijcnDaRxXCSjpy3id3NWU1JW4XU0CWGBVBxbgJbVtjOrHqtxjJlFAI2AgnpJJ+Khzs2SmD3pTK7t14rJ89dz5dOfsGGXlmgXbwRScSwGOphZGzOLAkYAs48YMxsYW/X1lcC7TqvESYiIjQrnweHdefq6PmwqOMSwv3zIjIWbtFCi1LuAKY6qcxaTgLnAamCWc26lmd1vZpdWDZsGpJhZLvBT4HtTdkWC3dBuzZh7xyD6tG7MPa+vYPz0HJ04l3plwf7bis/nczk5OV7HEKlzFRWO6Z9u5Hf//pKE6Aj+d3h3hnZr5nUsCRJmtsQ556vpuYDZ4xCR2gkLM24Y2IZ//fhMmjeK4ZaXlnDbK5+xR/c2Fz9TcYg0cB2aJvLGxIH8ZHBH5izfxpA/z+ftldu9jiVBTMUhEgQiw8O4fXAHZk86k/TEaCa8uIRJLy/VuQ/xCxWHSBDp2iKJNyYO5KdDOvL2yh0MfvgDXsvZrJlXUqdUHCJBJioijNvO78Cc28+iY9MEfvHXL7hu2kLW5x/wOpoECRWHSJBqn57AqxMG8MBlp/HF5n0MfeRDHn57je5xLqdMxSESxMLCjNEDsvjPz87mou7NePTdXC58ZD7vrdnpdTRpwFQcIiEgPSmGv4zozYzx/QgPM254bjE3Pr9Yy5bISVFxiISQge1Teev2QfxqWGcWbdjNBX/+gN/NWU1hUanX0aQBUXGIhJioiDAmDGrHuz8/m8t7ZTB5/nrOeeh9Xvx0I6XlWnVXjk/FIRKi0hNjeOiqnsyeNJD26Qn85s2VXPhI5cWDmr4rx6LiEAlxPTKTmTmhP8+MqVyWaMKLS7jq6U9ZuF53LJCaqThEBDNjSNemzL1jEA8O78bmPYe4ZsoCxj67iOV5+7yOJwFGq+OKyPcUlZbzwqcbefL9r9h7qJQhXZty+/kd6JbRyOtoUk+OtTquikNEjmp/USnPfbSRaR+tZ39RGYO7pHPb+R3okZnsdTTxMxWHikPklOwvKmX6xxuZ+tEG9h0u5awOqfzo7HYMaJeCmXkdT/xAxaHiEKkThUWlzFj4NdM+2kB+YTE9Mxtx89ntuKBrUyLCdco0mKg4VBwidaqotJy/L93C5PlfsangEJmNY7n+jCyuzm5JUkyk1/GkDqg4VBwiflFe4Zi3agfPfryBRRt2Ex8VzpV9Mrmuf2s6NE30Op6cAhWHikPE75bn7eO5jzfwzy+2UVJeQf+2TRjdP4shXZsSFaHDWA2NikPFIVJvCg4UMysnjxkLN5G35zAp8VEM753BNdkttRfSgKg4VBwi9a68wjF/bT6vLt7MO6t3UFbh6NUymStOz+DiHi1oEh/ldUQ5BhWHikPEU7sOFPPGZ1t4LSePNTsKiQgzzumUzuW9W3Be53TioiK8jihHCPjiMLMmwKtAFrARuNo5t6eGcW8B/YGPnHOXnMhrqzhEAodzjtXbCnlj2RbeXLaFHfuLiY0M57wu6VzSvTnndEonNirc65hCwyiOPwK7nXO/N7O7gMbOuTtrGHc+EAfcrOIQadjKKxwLNxQwZ/k23lqxnV0HSoiJDOOsDmkM6dKU87qkk5oQ7XXMkNUQimMNcI5zbpuZNQfed851OsrYc4CfqzhEgkdZeQWLNuxm7srtzFu1g637ijCDXi2TObtjGoM6ptEzM5nwMF2lXl8aQnHsdc4lV31twJ5vtmsYew7HKQ4zmwBMAGjVqlWfTZs21XFiEfEX5xyrtu1n3qodvL8mn8/z9uIcJMdFMqBtCme0S2FAuxTapSVouRM/CojiMLN3gGY1PHUPML16UZjZHudc46O8zjloj0MkZOw5WMKHubuYvzafT3J3sXVfEQBpidFkZzWmT+sm+Fo3pmuLJCK17EmdOVZx1NtUBufc4KM9Z2Y7zKx5tUNVO+srl4gEtsbxUVzaswWX9myBc46vdx/i068K+HR9ATkb9zBn+XYAoiPC6JbRiB6ZjeiZmUy3jEa0SY3X4S0/CJQ5cLOBscDvqz6/6W0cEQlEZkbrlHhap8Qzom8rALbvKyJn026WbtrLF3l7eWXR1zz38UYAYiLD6NQ0kS7Nk+jYNJH26Ql0aJpAs6QYHeY6BYFyjiMFmAW0AjZROR13t5n5gFucc+Orxn0IdAYSgAJgnHNu7rFeW4eqREJLWXkF63YeYMWWfazeVsjqbftZvX0/ew+VfjsmITqC1ilxZKXE0yoljtZN4shoHEuL5FhaNIptkFOCD5eUs7OwiO37ithRWMzO/UV0aJrI2R3TTur1AuIch1dUHCLinCP/QDG5Ow/w1c4D5O48wMaCQ3y9+xCbdx+irOK7/w42joskPTGGtMRo0hOjSU2MpnFcFI3jImkcH0Wj2EgSYyJIjK78HBsVTnRE2CnvxZRXOA6VlHGopJyDxWUcLC6nsKiUfYdL2V/1ec+hUvYcLGF31ceuA8XkFxZzsKT8e683sm8rfvfD7ieVJSDOcYiIeMXMSE+MIT0xhjPapX7nubLyCrbtK2Lr3sNs3XeYrXsrv95ZWPkP8oZdB8k/UExJWcVxfgbERoYTExlORJgRGR5GRLgRHmbYt2MM5xzlFY6yCkdZuaO0vILisgqKy8opLT/+L/KR4UbjuCiaxEfROC6K7pnJpCZEkZoQTVpiNM2SYmjWKIamSTEkxfjnn3gVh4iEtIjwMFo2iaNlk7ijjnHOcbi0nN0HS9hzsPK3/8KiMgqrPh8uLaeo6uNwaXlVITjKKir+uzdTrRO+KZRvCiYmsnKPJToinLiocOKjI4iPDicuKoKkmAiSYiMrP2IiSIiO8Pz8jIpDROQ4zIy4qAjioiLIrPFCgdCiSc8iIlIrKg4REakVFYeIiNSKikNERGpFxSEiIrWi4hARkVpRcYiISK2oOEREpFaCfq0qM8uncuHEk5EK7KrDOA2B3nNo0HsODafynls752pcITHoi+NUmFnO0Rb5ClZ6z6FB7zk0+Os961CViIjUiopDRERqRcVxbFO8DuABvefQoPccGvzynnWOQ0REakV7HCIiUisqDhERqRUVx1GY2VAzW2NmuWZ2l9d5/M3MnjWznWa2wuss9cXMWprZe2a2ysxWmtntXmfyNzOLMbNFZvZ51Xu+z+tM9cHMws3sMzP7p9dZ6oOZbTSz5Wa2zMxy6vz1dY7j+8wsHFgLDAHygMXASOfcKk+D+ZGZDQIOAC8457p5nac+mFlzoLlzbqmZJQJLgMuD/L+zAfHOuQNmFgl8BNzunFvgcTS/MrOfAj4gyTl3idd5/M3MNgI+55xfLnjUHkfN+gK5zrn1zrkSYCZwmceZ/Mo5Nx/Y7XWO+uSc2+acW1r1dSGwGsjwNpV/uUoHqjYjqz6C+rdHM8sELgamep0lWKg4apYBbK62nUeQ/4MS6swsC+gNLPQ4it9VHbZZBuwE5jnngv09PwL8EqjwOEd9csDbZrbEzCbU9YurOCTkmVkC8DfgDufcfq/z+Jtzrtw51wvIBPqaWdAemjSzS4CdzrklXmepZ2c6504HLgImVh2KrjMqjpptAVpW286sekyCTNVx/r8BM5xzf/c6T31yzu0F3gOGehzFnwYCl1Yd858JnGdmL3kbyf+cc1uqPu8EXqfy8HudUXHUbDHQwczamFkUMAKY7XEmqWNVJ4qnAaudcw97nac+mFmamSVXfR1L5QSQLz0N5UfOubudc5nOuSwq/z9+1zl3ncex/MrM4qsme2Bm8cAFQJ3OllRx1MA5VwZMAuZSecJ0lnNupbep/MvMXgE+BTqZWZ6ZjfM6Uz0YCIym8rfQZVUfw7wO5WfNgffM7Asqf0Ga55wLiSmqIaQp8JGZfQ4sAv7lnHurLn+ApuOKiEitaI9DRERqRcUhIiK1ouIQEZFaUXGIiEitqDhERKRWVBwiIlIrKg4REakVFYdIPTCzTDO75ojHssxsTtV9X9aa2d1e5ROpDRWHSP04Hzj9mw0zC6NyjaynnXOdgO6Azx8rmYrUNV05LuJnZnYm8CawFygEfgh0AsY7566oNq458IFzrqMXOUVOlPY4RPzMOfcRletCXeac6+WcWw90AT4/Ytw2IKlqYU2RgKXiEKkfnfjuKrTlQEL1AVWr9cYBZfWYS6TWVBwifmZmqcC+qlWXv/E+MKyqLL4xBFjqnAulO9VJA6TiEPG/LGBr9Qecc58DnwH3A5hZU+Bh4Ff1HU6ktlQcIv73JZBqZivM7AwAM7sL8AG/NrPzgKeA1sCTVfc/FwlYmlUlIiK1oj0OERGpFRWHiIjUiopDRERqRcUhIiK1ouIQEZFaUXGIiEitqDhERKRW/j+IA2InvlaT8AAAAABJRU5ErkJggg==\n", "text/plain": [ "
" ] }, "metadata": { "needs_background": "light" }, "output_type": "display_data" } ], "source": [ "t_1, z_1 = dynamics_1.expectations(0.5*sigma_z, real=True)\n", "plt.plot(t_1, z_1, label=r'$\\alpha=0.3$')\n", "plt.xlabel(r'$t\\,\\Omega$')\n", "plt.ylabel(r'$\\langle\\sigma_z\\rangle$')\n", "plt.legend()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Yay! This looks like the plot in figure 2a [Strathearn2018]." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Let's have a look at the above warning. It said:\n", "\n", "```\n", "WARNING: Estimating parameters for TEMPO computation. No guarantee that resulting TEMPO computation converges towards the correct dynamics! Please refer to the TEMPO documentation and check convergence by varying the parameters for TEMPO manually.\n", "```\n", "We got this message because we didn't tell the package what parameters to use for the TEMPO computation, but instead only specified a `tolerance`. The package tries it's best by implicitly calling the function `oqupy.guess_tempo_parameters()` to find parameters that are appropriate for the spectral density and system objects given." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "#### TEMPO Parameters" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "There are **three key parameters** to a TEMPO computation:\n", "\n", "* `dt` - Length of a time step $\\delta t$ - It should be small enough such that a trotterisation between the system Hamiltonian and the environment it valid, and the environment auto-correlation function is reasonably well sampled.\n", " \n", "* `tcut` (or `dkmax`) - Memory cut-off time (or number of steps). It must be large enough to capture all non-Markovian effects of the environment.\n", "\n", "* `epsrel` - The maximal relative error $\\epsilon_\\mathrm{rel}$ in the singular value truncation - It must be small enough such that the numerical compression (using tensor network algorithms) does not truncate relevant correlations." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "To choose the right set of initial parameters, we recommend to first use the `oqupy.guess_tempo_parameters()` function and then check with the helper function `oqupy.helpers.plot_correlations_with_parameters()` whether it satisfies the above requirements:" ] }, { "cell_type": "code", "execution_count": 11, "metadata": {}, "outputs": [ { "name": "stderr", "output_type": "stream", "text": [ "../oqupy/tempo.py:881: UserWarning: Estimating parameters for TEMPO computation. No guarantee that resulting TEMPO computation converges towards the correct dynamics! Please refer to the TEMPO documentation and check convergence by varying the parameters for TEMPO manually.\n", " warnings.warn(GUESS_WARNING_MSG, UserWarning)\n", "WARNING: Estimating parameters for TEMPO computation. No guarantee that resulting TEMPO computation converges towards the correct dynamics! Please refer to the TEMPO documentation and check convergence by varying the parameters for TEMPO manually.\n" ] }, { "name": "stdout", "output_type": "stream", "text": [ "----------------------------------------------\n", "TempoParameters object: Roughly estimated parameters\n", " Estimated with 'guess_tempo_parameters()'\n", " dt = 0.0625 \n", " tcut [dkmax] = 2.25 [36] \n", " epsrel = 2.4846963223857106e-05 \n", " add_correlation_time = None \n", "\n" ] } ], "source": [ "parameters = oqupy.guess_tempo_parameters(system=system,\n", " bath=bath,\n", " start_time=0.0,\n", " end_time=5.0,\n", " tolerance=0.01)\n", "print(parameters)" ] }, { "cell_type": "code", "execution_count": 12, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "" ] }, "execution_count": 12, "metadata": {}, "output_type": "execute_result" }, { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAX0AAAD4CAYAAAAAczaOAAAAOXRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjMuNCwgaHR0cHM6Ly9tYXRwbG90bGliLm9yZy8QVMy6AAAACXBIWXMAAAsTAAALEwEAmpwYAAAnRUlEQVR4nO3deZxU1Z338c+vqnqnu6GhodmkQQURZFFAMEaNBiTJBDRGMYuaJzEmJiZOYpKJM/NEEidPJotjNJlJNGqiEyeR7Ji4gMbEMWIQFRAEwQWh2bpZbJZeq+o8f9zqlWq6mq6t637fvupVVadO3fuzXs23b58691xzziEiIv4QyHQBIiKSPgp9EREfUeiLiPiIQl9ExEcU+iIiPhLKdAGdDRs2zFVXV2e6DBGRAeWFF17Y55yrTKRvVoV+dXU1a9asyXQZIiIDipm9lWhfDe+IiPiIQl9ExEcU+iIiPpJVY/r9sae+iarywkyXISJp1NraSk1NDU1NTZkuJS0KCwsZM2YMeXl5J7yNnAj91W8eYMldq1j26XnMrq7IdDkikiY1NTWUlpZSXV2NmWW6nJRyzrF//35qamoYP378CW8nKcM7ZnafmdWa2YZObUvNbKeZrY3d3puMfXUXjkS56VfrcMBNy9YRjkRTsRsRyUJNTU0MHTo05wMfwMwYOnRov/+qSdaY/s+AhXHab3fOzYjdHknSvrq4/9lt7DvSDMC+I808sGpbKnYjIlnKD4HfJhn/r0kJfefc08CBZGyrL2oPN3Hbyi00tkQAaGiJ8L0VW6g73JzuUkREBoRUz965wczWx4Z/hsTrYGbXmdkaM1tTV1fXp40vX7uLSLTr9QAiUcfydbtOvGIRkTSqrq5m3759adtfKkP/R8DJwAxgN3BbvE7Oubudc7Occ7MqKxM6i7jd4hmjCQa6/rkTMFg0fdQJFSwi0h/OOaLR7P5eMWWh75zb65yLOOeiwE+AOcneR2VpATfNn0hRfrC9beGUKipLC5K9KxGRuLZt28akSZO4+uqrmTp1KrfeeiuzZ89m2rRp3HLLLe39LrnkEs466yymTJnC3XffnbF6UzZl08xGOud2x55eCmw4Xv8Tdc051dy/6i22H2ggL2is31mPc85XX+6ICHz94Y28sutQUrd5+qgybnn/lF77bd26lfvvv59Dhw7x61//mtWrV+OcY9GiRTz99NOcd9553HfffVRUVNDY2Mjs2bO57LLLGDp0aFLrTUSypmz+AlgFTDKzGjP7BPAdM3vZzNYD7wK+kIx9dRcKBvje5dMxvF8Ar9cdpeZgYyp2JSIS17hx45g7dy4rVqxgxYoVzJw5kzPPPJPNmzezdetWAO68806mT5/O3Llz2bFjR3t7uiXlSN8596E4zfcmY9uJmDO+glU3X8Su+kbu+d832bznMGMritO1exHJAokckadKSUkJ4I3p33zzzXzqU5/q8vpf/vIXnnjiCVatWkVxcTEXXHBBxs4izpm1d6rKC5k4ohSAV/ck9088EZFEXHzxxdx3330cOXIEgJ07d1JbW0t9fT1DhgyhuLiYzZs389xzz2WsxpxYhqHNoIIQYyuK2LzncKZLEREfWrBgAZs2bWLevHkADBo0iJ///OcsXLiQH//4x0yePJlJkyYxd+7cjNWYU6EPMGlEKa8q9EUkTaqrq9mwoWOeyo033siNN954TL9HH3007vu3bduWqtLiypnhnTaTqkp5Y99RmsORTJciIpJ1cjD0y4hEHW/UHc10KSIiWSfnQv+0qrYvczXEIyLSXc6F/vhhJeQFTV/miojEkXOhnxcMcHLlIE3bFBGJI+dCH7wvczW8IyJyrJwM/eqhJew+1ERLOLtXuxORge+cc87JdAl9kpOhP3pIEc55F0sXEUmlZ599NtMl9ElOhv6YwUUA1LzdkOFKRCTXDRo0CPDW1zn//PNZvHgxEyZM4Ktf/SoPPvggc+bM4YwzzuD1118H4OGHH+bss89m5syZvPvd72bv3r0A1NXVMX/+fKZMmcK1117LuHHjUnJxlZw7Ixe8I32AnVptU8Q/Hv0q7Hk5udusOgPe8+8Jd1+3bh2bNm2ioqKCCRMmcO2117J69WruuOMOfvCDH/D973+fc889l+eeew4z45577uE73/kOt912G1//+te58MILufnmm3nssce4997UrFmZk6FfVV4IwM63Ffoikj6zZ89m5MiRAJx88sksWLAAgDPOOIOnnnoKgJqaGpYsWcLu3btpaWlh/PjxADzzzDP87ne/A2DhwoUMGRL3CrP9lpOhXxAKMry0QEf6In7ShyPyVCko6LhqXyAQaH8eCAQIh8MAfO5zn+OLX/wiixYt4i9/+QtLly5Na405OaYP3hCPjvRFJNvU19czevRoAO6///729ne84x0sW7YMgBUrVnDw4MGU7D93Q3+wQl9Ess/SpUu5/PLLOeussxg2bFh7+y233MKKFSuYOnUqv/rVr6iqqqK0tDTp+zfnXNI3eqJmzZrl1qxZk5RtfevRTfz0mW1svnUhgYCulyuSizZt2sTkyZMzXUZSNDc3EwwGCYVCrFq1iuuvv561a9ce0y/e/7OZveCcm5XIfnJyTB+8aZstkSh1R5oZUVaY6XJERI5r+/btXHHFFUSjUfLz8/nJT36Skv3kbOi3Tdtcv6Oe+VMU+iKS3U499VReeumllO8nZ8f0R8VO0Prkf6/h+W0HMlyNiKRKNg1Rp1oy/l9zNvRHlHZMnbpp2TrCEa3DI5JrCgsL2b9/vy+C3znH/v37KSzs38hFzg7v/PbFne2P9x1p5oFV2/j4uRMyWJGIJNuYMWOoqamhrq4u06WkRWFhIWPGjOnXNnIy9GsPN3Hbyi3tzxtaInxvxRbeP300lZ3+AhCRgS0vL6/9jFZJTE4O7yxfu4tItOufe5GoY/m6XRmqSEQkO+Rk6C+eMZpgt7n5wYCxaPqoDFUkIpIdcjL0K0sLuGn+REJBL/iLQgG+tGCihnZExPdyMvQBrjmnmrLCPAAGl+Rz9bzqzBYkIpIFcjb0Q8EA159/MgCfPv9kQsGc/V8VEUlYTifh+ZMqARhcnJfhSkREskNOh37bmjt7D+lauSIikKTQN7P7zKzWzDZ0aqsws5VmtjV2n5rLwBxHWWGIorwgew81p3vXIiJZKVlH+j8DFnZr+yrwpHPuVODJ2PO0MjOqygvZoyN9EREgSaHvnHsa6L6q2WKg7bIw9wOXJGNffTWirIC99Qp9ERFI7Zj+COfc7tjjPcCIeJ3M7DozW2Nma1KxfkZVmY70RUTapOWLXOctgRd3GTzn3N3OuVnOuVmVlZVJ3/eI8kJqDzX7YhU+EZHepDL095rZSIDYfW0K99WjqrJCWiJRDhxtycTuRUSySipDfzlwTezxNcAfUrivHlXFpm1qiEdEJHlTNn8BrAImmVmNmX0C+HdgvpltBd4de552wzVXX0SkXVLW03fOfaiHly5Kxvb7o6q8LfQ1V19EJKfPyAUYXlqAGezRtE0RkdwP/bxggKElBRreERHBB6EPUFVeoC9yRUTwS+iXFWp4R0QEn4T+iLJCDe+IiOCT0K8qK+RgQytNrZFMlyIiklG+CP0RsWmbtZq2KSI+54/Q11m5IiKAT0K/SmfliogACn0REV/xReiXFYUozAto2qaI+J4vQt/MdDEVERF8EvqgufoiIuCj0NcF0kVE/BT6ZYXs1WUTRcTnfBP6w8sKaQlHOdjQmulSREQyxjeh337ZRM3gEREf80/olxcAsPewQl9E/Ms3od+2FMNeHemLiI/5JvSHl2r9HRER34R+fijAsEEFGtMXEV/zTegDjBpcyC6Fvoj4mL9Cv7yIXW83ZroMEZGM8VXojxxcyO63G3WCloj4lq9Cf/TgIo62RDjUGM50KSIiGeGr0B9ZXgTArnoN8YiIP/kr9Ad70zY1ri8ifuWr0B89uO1IXzN4RMSffBX6wwYVEAqYjvRFxLd8FfrBgFFV7s3gERHxI1+FPrTN1dfwjoj4U8pD38y2mdnLZrbWzNaken+98c7KbdRyDCLiS+k60n+Xc26Gc25WmvbXo5GDi9hd38Tcbz3J89sOZLocEZG08t3wTlVZAZGod0buTcvWEY5EM1yRiEj6pCP0HbDCzF4ws+u6v2hm15nZGjNbU1dXl/JiXtl1qP3xviPNPLBqW8r3KSKSLdIR+uc6584E3gN81szO6/yic+5u59ws59ysysrKlBZSe7iJP6zd1f68oSXC91Zsoe5wc0r3KyKSLVIe+s65nbH7WuB3wJxU77Mny9fuItptsbVI1LF83a4e3iEikltSGvpmVmJmpW2PgQXAhlTu83gWzxhNKNj1fzkYMBZNH5WhikRE0ivVR/ojgGfMbB2wGviTc+6xFO+zR5WlBdw0fyJm3vPivCBfWjCRytKCTJUkIpJWoVRu3Dn3BjA9lfvoq2vOqea2lVtoaIkwrLSAq+dVZ7okEZG08d2UzVAwwILTqwD4zgenHTPcIyKSy3yZeHMnVAAdq26KiPiFL0P/pKHFAGw/0JDhSkRE0suXoT9uaAkAb+1X6IuIv/gy9KvKCskPBnjrwNFMlyIikla+DP1gwBgzpIjtOtIXEZ/xZeiDN66vMX0R8Rvfhv64imK272/AdVuWQUQkl/k29MdWFHO4OczBhtZMlyIikja+Df22GTwa4hERP/Ft6FfH5uq/ue9IhisREUkf34b+uKElhALG1r0KfRHxD9+Gfn4owPhhJWxR6IuIj/g29AEmjihla+3hTJchIpI2vg79U0cMYvuBBhpbIpkuRUQkLXwd+hNHlOIcvF6nIR4R8Qefh/4gALbs1RCPiPiDr0N/3NAS8oKmL3NFxDd8Hfp5wQAThg1iq470RcQnfB364H2Zu0UzeETEJ3wf+hNHlLLjQCMNLeFMlyIiknIK/RGlABrXFxFf8H3oTxtTDsC6HW9nthARkTTwfeiPLC9kRFkBL20/mOlSRERSzvehb2bMHDuEl3SkLyI+4PvQB5hx0mDe2t/A5t2HMl2KiEhKKfSBmWMHA7Dwjv/l+W0HMluMiEgKKfSBySNL2x/ftGwd4Ug0g9WIiKSOQh/41ZoazLzH+44088CqbRmtR0QkVXwf+rWHm7ht5Rac8543tET43oot1B1uzmxhIiIp4PvQX752F5Go69IWiTqWr9uVoYpEPHvqmxJq62t7T33FH1Ie+ma20MxeNbPXzOyrSd/B+mVw+1RYOti7X7+sT29fPGM0wYB1aQsGjEXTRyWxSMm0dAdof7ex+s0DzPvWk10mFsRr67F9/TKav3s6w/9jBM3fPd37dxGvTXzHnHO99zrRjZsFgS3AfKAGeB74kHPulXj9Z82a5dasWZP4DtYvg4c/D62NHW15RfD+O2HaFQlv5t7/fYPvrdxCY0uEgMG/vm8yHz93QuJ1SOLWL4Mnv4Grr8HKx8BFX/Pau7dNu6L/fWPtzY8vJe/ITloHjabg4qUAx7b11DcD24isfYi6P/wLw6N11AYqqVz8TYBj2oIzlsTtGwwY7uHPY53+XbhAHphhkZaOtrwirI//ViQ7mdkLzrlZCfVNcejPA5Y65y6OPb8ZwDn3rXj9+xz6t0+F+h3HtpePhS9sSHgz4UiUC2/7K9sPNACw+p8vYnhZYeJ15JpYgFJfA3GCtUs79C2wu/+SjoURncKIvCKY/mFY9z8n3jfW7tb9T2LhF69vqBA7Ywnu5YewcMfRuLOQt41oa9e+kxfjNv0eCzd3ai/ATvsH3OY/dm23EBhYNNy176kXE978GCHXUV+YIAYE6bisZ9jyCU2cT3jLyq59LY9QXiG0JLhybPEwuORHbRXE7ly3x51ei9d23D4JtPeUQT1to0+vdXmhl3300rfH/r33dd3e6nDeu1zHcxsxldBJs3uo5/iyKfQ/CCx0zl0be34VcLZz7oZOfa4DrgM46aSTznrrrbcS38HSwcT/wA2Wvt2nWle/eYAr7loFwPeXzOCSmaP79P6sd5zAPqYt3l9PiYZwqBCmfAA2/gY6hRzBfO/WkujCdkaP/5jkhO2IDmNsYF+my5A4Vo28mnmf+sEJvXdAhX5nmTrSb7Pr7UYW/fAZzjl5GHd+aGaf358Vegr3RI6wQ4UQLIDm+jgbzo4Qds6rIkqAKEaUAJHuNxckQoAwHffttzivtRIkTMh77oK0tj3u9L4IAVoJddtGiHCnbbUSwmHezVl7nd5xc9u993q0rR9A7L5re+fXvO+cOr+/7THd7p3r3pcufaLOMPP+ajCc96VeXiFmAcy8ZUkC3gPMDMNrMzMsYHj/tbUBFiAQex/QsR2sfXvevff1YcDr1NEntk3a+nbaH3S0dTw2LEDsPW19uu0PcG19vU239+l4W2xfnWqm/XnHk87b6KjR2l+jvWvH/1NbTW0POj1s70fntljD5HGjOHvy+J5/+I+jL6EfOqE9JG4nMLbT8zGxtuS46Gvxj0rbgq6PRg0u4vyJw1n5yh6aWiMU5gWTVGgKJBLu9Tu858HCrp8RQKehiXbhJu8W14kFfosL0kABDa6QIxRy1BVxGO92xBVziGKOuCIOUcxRijjqCmnAuzW6fBopoJF8mpx330gB0X7OPwgRJp8w+bSST5g8wuRZmDwi3mO8x/nWSpAo+bF4b3s9ZBFCRMgjQpAIhbQSJEKQKHkWJkSUkHmv5cXCNZ9WArhYe7S9fyh2HyAau3fkWbhbWzT2Hu/1tsdmxH4Neb/+grHXogWlRJqOUmwt7e+NOC98CgIdw0mNLp/W991B2ZwP9+vzlIEl1aH/PHCqmY3HC/srgeT9hLV9AbXyFji8CwoHw3u/268vpi6dOZrfvFjD4xv3sHhGlg7xdP8C+3jh3tp4bNsJCJPHQVfEflfOAVfKAUo52Om+3pXwNoOodyXUU0K9K+EQxbSQ3+u2A0QpoZFimimxJgppIRwqZkR4D6NsP8U0UxRoJp8WCmllkDVRQAtFNBPA8VbR6ZzR+Dxl1kCBtVJAKzioGX4ep9U9Rpk1kE841u4IGORb1/BrnLKEoo0PUWQtvba3uCBgKdmGm/ER8l7+BXnRjl++YUKYGUHX8Yu6NVBI3pkfofXFB7v0bQ0UkvcPt/LU5lombrydkexnN0PZesYXcI6ubVO/wAUKfN9Jaeg758JmdgPwOBAE7nPObUzqTqZdAVMuhVuHwdzP9HsmwjknD2Xc0GIe/Pv27An97l+KthxNWrhHnbGPcna5oexxFewOjmZXazH7XRl1DKbODabWDeYgpe1DBN2V0kC5HWEIRxhEA0PyCzmzdS0Vdpgya6CUBgpp4c0xi9hUs5+reYTxgd0ccUXcxyICRPmC/ZJRtp9dbih3cCWfuP6fuPdH3+ZGOtq/z5UcBj5Ct74fP7bvHXYlX/nY19jw6KlM3Hg7Qzl8/PC7/AaeWnZawu0p28alNxAZP489f/gXhkf3URsY1j57p0vbom/CjCUExpx9bPu0K3jnlCgXbpvC9gMNnFRRzJ8/cD5A3Dbxl1Qf6eOcewR4JKU7CeZBQTk09n+xtEDA+NCck/j3RzfzWu1hThle2vubUqnTUb1B/O8wenEwfxRbG0vZTQU73HC2u+Fsjw5npxvGHobS2vnHoBWMKMN5myo7QCVv0xzIp7xiJOUHXuZS+yunBHYRcQF+6S4ijzBfDnYL4U9fz70/+jZXsbJL+1c+chXL1+7ksyun09gSoTgvyJcunohzMH/lO2ls7mg7bWQZp83/xDHtfelbWVrAOy/7TMLh11PfdG8jOGMJ28vnM++uVSz79DyqqisA4rb11DcUDPC9y6ez5K5V3HbFdEJBb0gsXpv4S0q/yO2rPn+R29kdM2DMLLjsnn7Xse9IM/O+9SRXza3ma+8/nT31TVSVp2EKZ/cj+nf9C6z4v9BQ1+tbnYOd+ePZ2FjBdkbwmhvNa9FRvOFGcpCyLn2HUo/DKLZmLrC1nBqooYBWVnA2p77jUh54dhuNrR2LzhXnB/nt9efwgR89S0NLpEv7p86bwI+ffqNLiH/83Aldzn3o3N55euxJFcX8+aZY+HVrCwUD/e7bFmqr3zzAklgozo6FYry2vranchtA3J+7nn4W+9Ketp9nSZusmb3TV/0K/Z9cCIXlcNXvklLL53/xEk9s2svtV8zg0z9/Ie4/ymRo/wcY70SzHkScsTU6mg2MZ0N0PK9Ex7HZncQhStr7DKWesVbL25Qy9pSp/P3NA7SEO4K8KC8AGI2tXUP8r19+F8vX7owb2PGC/Op51WkP4b6GbZfPuZe2vrYrhCUb+DP0H7wcjtTCp/6alFq27TvK/Nv/Sl4wQENL5JjgSoYuIfXbd8KhmmP6OALscyW8ED2Nl6Kn8FL0FF5242nEC4kimjjZdnPUihg/cRp/e20fzeGuR+k9BblzJHw0frwgz0QI9zVsRXJZX0I/dwb1iodCQ/IugFI9rIRZ4yrahzP6u+Ry9/VVImsf4qQH5vB6wYc5+WczcJ0C/7Ar4pHIbP6p9ZO8q/m7zG6+i0+3foGfRhbSQojLgk9zVmALAI0UssGNZ29oNN++bBpfuXgSRfneVNPivCBfWuCNbV9zTjWVgwoAGFZawNXzquO2Qcd4sEGXsd+e2ueMr2DVzRcdE+w9tQNxg7mnsO5vXxHpkDuhX1SRlC9y29QebmLtjo6LpcdbcvmEF89av4zo8s9T5eoIGFRQT5PL5/fhc1jQ/G3OaL6Xz7R+gcciszk57wDTgtv4Ueh21uVfyw/zfkBzcBDnXfi+hMMd4gd2TyEOfQ9yhbDIwJA7oV9c4Z3i35qcZWOXr91FtxWXCUei7UsuJ7oKYjgS5aZfrcPRcVWu1se/zsFIARFnHHUF/Fvrh5ne/BP+OXwtY6yOfw39nD/l38yz+Z/jo2eUsHjhe/hi4EtMbrmfBe4/mbzgE3z2Xaf06Sgd4gd2X4/Gj9cuItkvh0J/qHefpKP9eEsuR6KOCyZWxg3yeG0A9z+7jX1HvL8O9tQ38o8/+CVvHnYMt3r+EH0H17V8kXwiPJD/LV7M/yTfCP2MjwcfpZyjLOWTTF14bdyj9xM5Soe+DZWISO5J+Tz9tCmOBVzDASjr/1r4laUF3DR/YvsXnXkBI+IcV9+3mjnVQ6g77P1F0TbW7xzt4b7vSDN3/fV1Rg4uYtmjTzDOhXmVsbzXPcv/O3gPEQvxt8jpnBt4mQ8UPNO+z0MFVcxvvOOYueYQf351W7gryEUkUbkze2fbM/Cz98HVy2FCcs407D5b5ftLpvPFZevYtr+hS79Q7C+CcPfxICCPMGfaFr6Ut4zZgS1EnPHayPcxofaJY06fDyy6k3etGB53miNoZoqIxOff2TsADfuTtsnuQyhnjqvgI2efRF63YZ9w1B0T+MGAcW/Jf7Ku4FoeKvg3Zsdm2wTNccrRtQQW3ckeqyTqjD1WSWDRnQRnLOlxyAZ09C4i/Zc7oV8UG95J4gweOHZ8/JKZY8gLdf3YCkNGYbe2glCACyN/o7jTIlttgod3eqfPX72ak5sfZMc1qwnOWBJ3fyIiyZQ7od95TD/JOh9ht431d54u+ZWFp/HlbvPjv/zuCVioKP4Gy8cAfZ/+KCLSX7kT+m2LrqUg9Lvr9USnQflcU/9fEG70LlbSWbf1/hXwIpJOuRP6AMVDkjqm35O4Jzpt/DUr7bO8UfBh/hy+hsALP4V3/CNc8l9QPtZblrh8bJ8v2i4ikky5M2UTYksxpD70odt0ydhiaQWtjWAQaD0EFoThk72An3ZFDyvRi4ikV44d6Q9N+he5x9M+NPPkN45dHdNF4M//lrZaREQSkVuhX1SRtiP9LuqPXR3zuO0iIhmSW6Gf5JU2ExabjZNwu4hIhuRY6McWXQs39943mc7+9LFt3WbpiIhkg9wK/ZJh3v3Rfenbp3PwxlMQKoTSkaBZOiKSxXJr9k5JpXd/tA7KR6dnn68+Cq89AQu+CefckJ59ioicoBw70h/u3R/t/ULiSdHaBI/fDMMmwdmfSs8+RUT6IceO9NuGd9IU+s/fAwe3wVW/984IFhHJcjkW+rHhnSO1qd3P+mXwxFI4tBNCBen7JSMi0k+5FfoFgyCvOLUhHDv7tv1krHCz9xz0xa2IZL3cGtMHb4gnlaEf7+zb1kavXUQky+Vg6A9Pbejr7FsRGcByMPQrUxv6ZT1MBdXZtyIyAORg6A+DIykM/VPnH9ums29FZIDIvdAfNBwa9kE0mvxtOwc7VntH+2Vj0Nm3IjLQ5NbsHfCGd6JhaHq74xKKyfLWs1C7ERb9EM68KrnbFhFJg5Qd6ZvZUjPbaWZrY7f3pmpfXXReiiHZXvipd0nGqZclf9siImmQ6uGd251zM2K3R1K8L0+qQv/ofnjlDzB9CeQXJ3fbIiJpkntj+qk6K3fdLyDSAmd9LLnbFRFJo1SH/g1mtt7M7jOzIfE6mNl1ZrbGzNbU1SXh6HxQ26JrSVxe2Tl44WcwZg6MmJK87YqIpFm/Qt/MnjCzDXFui4EfAScDM4DdwG3xtuGcu9s5N8s5N6uysrI/5XiKhoAF4GgSj/R3vQj7t8LMjyZvmyIiGdCv2TvOuXcn0s/MfgL8sT/7Slgg6F02MZlj+usegmABnL44edsUEcmAVM7eGdnp6aXAhlTt6xglw5Nzgtb6ZfAfU2D1Xd5fD1tX9H+bIiIZlMp5+t8xsxmAA7YB6bvKSDIWXTtmNc1GraYpIgNeykLfOZe5s5cGDYeaNf3bxvFW01Toi8gAlXtTNsEb3jm8x5t1c6K0mqaI5KDcDP3yMd5wTOPB/m2jL+0iIgNAjoZ+bPnj/hyVX/Q1sGDXNq2mKSIDXG6GflnsaPzQzhPfxumLvYud55Wg1TRFJFfk3iqb0DEE058j/defgnATfOQ3cGpCpyOIiGS93DzSL6mEQF7/Qv+VP0BhOYw/L3l1iYhkWG6GfiDgjeufaOhHwrDlUZi4EEL5ya1NRCSDcjP0wRvXP9Ex/R3PeTN/JqXnEgAiIumSu6FfPubEj/RffRSC+XDKRcmtSUQkw3I49EfDoV0QjfTtfc7B5j/B+POhoDQ1tYmIZEgOh/4YcBHvzNy+qNsMB9+E0zS0IyK5J3dD/0Tn6m/+k3c/8T3JrUdEJAvkbui3z9Xf0bf3bXkMRp0JZSN77ysiMsDkcOi3LcXQhyP9I7Xe6pyTdJQvIrkpd0N/y+OAwcr/C7dP9dbHT+g9zpufLyKSg3JzGYa2C6AQW1q5fsfxL4Cyfpm3Tn79Dm+RtdpNMHJa2soVEUmX3DzSP94FULpr+wXRNvbvIvDHGxP7y0BEZIDJzdDvywVQ+vILQkRkgMvN0O/LBVB0hSwR8ZHcDP2LvuZd8KSzni6AoitkiYiP5GboT7vCu+BJ2wlaBaU9XwDloq956+x0pitkiUiOys3QBy/gv7gRhk+Bk+b1fMWraVd46+y00RWyRCSH5eaUzc5GTIG3/tbz687Bvi1wynz46K/TV5eISAbk7pF+m6qp3vo7DQfiv753A7z9Fkz+h/TWJSKSAbkf+iOmePd7N8Z/fd0vIRCC0xT6IpL7fBD6Z3j38UI/0grrH/KWXSgZlt66REQyIPdDf9BwKB4Ge9Yf+9rWFXC0DmZ+NP11iYhkQO6HvhmMPw9efQRam7q+9tKDMGiE9yWuiIgP5H7oA5x5tXeh800Pd7TV7/TWzp+2BIK5P4lJRAT8MGUTvHn4Q6rhxfsB17GiJgalVRkuTkQkffwR+oGAd7T/5DegZjWEm2MvOPjzrVBSqZOxRMQX+jW8Y2aXm9lGM4ua2axur91sZq+Z2atmdnH/ykyCGR/x7tsDP0YraoqIj/T3SH8D8AHgrs6NZnY6cCUwBRgFPGFmE51zkX7u78QdbxhHK2qKiE/060jfObfJOfdqnJcWA790zjU7594EXgPm9GdfSVE+tod2ragpIv6Qqtk7o4EdnZ7XxNqOYWbXmdkaM1tTV1eXonJi+rLksohIDuo19M3sCTPbEOe2OBkFOOfuds7Ncs7NqqysTMYme9a25HL5WMC0oqaI+E6vY/rOuXefwHZ3Ap3HUsbE2jJv2hUKeRHxrVQN7ywHrjSzAjMbD5wKrE7RvkREJEH9nbJ5qZnVAPOAP5nZ4wDOuY3AMuAV4DHgsxmduSMiIkA/p2w6534H/K6H174JfLM/2xcRkeTyx9o7IiICKPRFRHzFnHOZrqGdmdUBb53g24cB+5JYTroMxLpVc3qo5vTIhZrHOecSmvOeVaHfH2a2xjk3q/ee2WUg1q2a00M1p4ffatbwjoiIjyj0RUR8JJdC/+5MF3CCBmLdqjk9VHN6+KrmnBnTFxGR3uXSkb6IiPRCoS8i4iMDLvTNbGHsEoyvmdlX47xeYGYPxV7/u5lVZ6DM7jX1VvPHzKzOzNbGbtdmos5uNd1nZrVmtqGH183M7oz9P603szPTXWOcmnqr+QIzq+/0OWf8QgpmNtbMnjKzV2KXHr0xTp+s+qwTrDmrPmszKzSz1Wa2Llbz1+P0yarsSLDmvmeHc27A3IAg8DowAcgH1gGnd+vzGeDHscdXAg8NgJo/Bvww059vt5rOA84ENvTw+nuBRwED5gJ/HwA1XwD8MdN1dqtpJHBm7HEpsCXOz0dWfdYJ1pxVn3XssxsUe5wH/B2Y261PtmVHIjX3OTsG2pH+HOA159wbzrkW4Jd4l2bsbDFwf+zxr4GLzMzSWGN3idScdZxzTwMHjtNlMfCA8zwHDDazkempLr4Eas46zrndzrkXY48PA5s49ipzWfVZJ1hzVol9dkdiT/Nit+6zWLIqOxKsuc8GWugnchnG9j7OuTBQDwxNS3XxJXrpyMtif7r/2sx6uJhvVkn4kphZZl7sz+VHzWxKpovpLDacMBPviK6zrP2sj1MzZNlnbWZBM1sL1AIrnXM9fs5Zkh2J1Ax9zI6BFvq56mGg2jk3DVhJx9GGJNeLeGuUTAd+APw+s+V0MLNBwG+Af3TOHcp0PYnopeas+6ydcxHn3Ay8K/nNMbOpGS6pVwnU3OfsGGihn8hlGNv7mFkIKAf2p6W6+Hqt2Tm33znXHHt6D3BWmmrrj+y9JGYPnHOH2v5cds49AuSZ2bAMl4WZ5eGF54POud/G6ZJ1n3VvNWfrZw3gnHsbeApY2O2lbMuOdj3VfCLZMdBC/3ngVDMbb2b5eF+2LO/WZzlwTezxB4E/u9g3HhnSa83dxmcX4Y2RZrvlwNWxmSVzgXrn3O5MF3U8ZlbVNkZrZnPwfv4z+o86Vs+9wCbn3H/00C2rPutEas62z9rMKs1scOxxETAf2NytW1ZlRyI1n0h29OvKWenmnAub2Q3A43izYu5zzm00s28Aa5xzy/F+GP/bzF7D+1LvysxVnHDNnzezRUAYr+aPZazgGDP7Bd4MjGHmXRLzFrwvknDO/Rh4BG9WyWtAA/B/MlNphwRq/iBwvZmFgUbgygwfEAC8A7gKeDk2dgvwz8BJkLWfdSI1Z9tnPRK438yCeL+Aljnn/pjN2UFiNfc5O7QMg4iIjwy04R0REekHhb6IiI8o9EVEfEShLyLiIwp9EREfUeiLiPiIQl9ExEf+PwiQXzLi2o4AAAAAAElFTkSuQmCC\n", "text/plain": [ "
" ] }, "metadata": { "needs_background": "light" }, "output_type": "display_data" } ], "source": [ "fig, ax = plt.subplots(1,1)\n", "oqupy.helpers.plot_correlations_with_parameters(bath.correlations, parameters, ax=ax)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "In this plot you see the real and imaginary part of the environments auto-correlation as a function of the delay time $\\tau$ and the sampling of it corresponding the the chosen parameters. Sampling points with spacing `dt` are given up to time `tcut` (`dkmax` points). We can see that the auto-correlation function is close to zero for delay times larger than approx $2 \\Omega^{-1}$ and that the sampling points follow the curve reasonably well. Thus this is a reasonable set of parameters." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We can choose a set of parameters by hand and bundle them into a `TempoParameters` object," ] }, { "cell_type": "code", "execution_count": 13, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "----------------------------------------------\n", "TempoParameters object: my rough parameters\n", " __no_description__\n", " dt = 0.1 \n", " tcut [dkmax] = 3.0 [30] \n", " epsrel = 0.0001 \n", " add_correlation_time = None \n", "\n" ] } ], "source": [ "tempo_parameters = oqupy.TempoParameters(dt=0.1, tcut=3.0, epsrel=10**(-4), name=\"my rough parameters\")\n", "print(tempo_parameters)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "and check again with the helper function:" ] }, { "cell_type": "code", "execution_count": 14, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "" ] }, "execution_count": 14, "metadata": {}, "output_type": "execute_result" }, { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAXkAAAD4CAYAAAAJmJb0AAAAOXRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjMuNCwgaHR0cHM6Ly9tYXRwbG90bGliLm9yZy8QVMy6AAAACXBIWXMAAAsTAAALEwEAmpwYAAAj6UlEQVR4nO3deZxU1Z338c+vqrqrVxp6UXYbjSACgtogomOMCyIz0SzjlkXnGR2dbDJGk0cnz6hxJtGgxmgmG0YmJJqFJOMTSaLBGB2jYhATQCIILgjNIt00tkDvVWf+uNUbXd30Uktz7/f9evWrqu65y6+L7m9fTp17rjnnEBERfwpluwAREUkfhbyIiI8p5EVEfEwhLyLiYwp5EREfi2S7gK7Ky8tdZWVltssQETmivPzyy7XOuYpkbcMq5CsrK1mzZk22yxAROaKY2du9tam7RkTExxTyIiI+ppAXEfGxYdUnPxS765sYXZKX7TJEJINaW1uprq6mqakp26VkRF5eHuPHjycnJ6ff2/gi5Fe/Vcdl31vF8n8+ndmVpdkuR0QypLq6muLiYiorKzGzbJeTVs459u7dS3V1NZMmTer3dinprjGzpWa2x8w2dFl2u5ntMLO1ia+FqTjWodpicW78+ToccOPydbTF4uk4jIgMQ01NTZSVlfk+4AHMjLKysgH/ryVVffI/ABYkWX6fc25W4uu3KTpWN8te2ErtgWYAag8088NVW9NxGBEZpoIQ8O0G872mJOSdc88CdanY10Ds2d/EvU9uprElBkBDS4x7Vm6mZn9zpksRERmW0j265rNmtj7RnTMq2Qpmdq2ZrTGzNTU1NQPa+WNrdxKLd58PPxZ3PLZu5+ArFhHJoMrKSmpra9O2/3SG/HeA44BZwC7g3mQrOeeWOOeqnHNVFRVJr8rt1cWzxhEOdf/vSzhkXDRz7KAKFhEZCucc8fjw+lwwbSHvnHvHORdzzsWBB4E5qT5GRXGUG8+fTH5uGICQwU3zJ1NRHE31oUREktq6dStTpkzhyiuvZPr06fz7v/87s2fP5qSTTuK2227rWO9DH/oQp556KtOmTWPJkiUZqy9tQyjNbIxzblfi5YeBDX2tP1hXzatk2aq32VbXgAM+ftox6TiMiAxzX17xV17d+V5K93ni2BHc9sFph11vy5YtLFu2jPfee49f/OIXrF69GuccF110Ec8++yxnnXUWS5cupbS0lMbGRmbPns1HP/pRysrKUlpvMqkaQvkTYBUwxcyqzexqYLGZvWJm64EPADek4liHioRD3HPJTACcg9drDqTjMCIivTrmmGOYO3cuK1euZOXKlZx88smccsopbNq0iS1btgDwwAMPMHPmTObOncv27ds7lqdbSs7knXNXJFn8UCr23R9zJpXy6Kfn8eFvv8CarfuYNrYkU4cWkWGiP2fc6VJYWAh4ffK33HIL1113Xbf2Z555ht///vesWrWKgoICzj777IxdpeubuWtmTRjJmJI8Xtqa8ZGcIiIAXHDBBSxdupQDB7wehR07drBnzx7q6+sZNWoUBQUFbNq0iRdffDFjNfliWgPwLhKoqixl9Vt7cc4F6gIJERke5s+fz8aNGzn99NMBKCoq4uGHH2bBggV897vfZerUqUyZMoW5c+dmrCbfhDzAnMpRrFi3k+p9jUwoLch2OSISAJWVlWzY0DmuZNGiRSxatKjHeo8//njS7bdu3Zqu0gAfddcAVCUmJ1vztrpsRETAZyF/bIX34Ud1XWOWKxERGR58FfLRSJiS/BxqDmjuGhER8FnIg3cVrCYoExHx+C/kixTyIiLt/BfyxVF114iIJPgz5HUmLyIZMm/evGyX0CdfhnxDS4yDzW3ZLkVEAuCFF17Idgl98l/IF3nTDNeqy0ZEMqCoqAjw5qd5//vfz8UXX8yxxx7LzTffzCOPPMKcOXOYMWMGb7zxBgArVqzgtNNO4+STT+a8887jnXfeAaCmpobzzz+fadOmcc0113DMMcek5GYivrriFeiYS75mfzPHlBVmuRoRyZjHb4bdr6R2n6NnwIV39Xv1devWsXHjRkpLSzn22GO55pprWL16Nffffz/f/OY3+cY3vsGZZ57Jiy++iJnx/e9/n8WLF3Pvvffy5S9/mXPOOYdbbrmFJ554goceSs0cj74OeRGRTJo9ezZjxowB4LjjjmP+/PkAzJgxg6effhqA6upqLrvsMnbt2kVLSwuTJk0C4LnnnuPRRx8FYMGCBYwalfSOqQPmu5AvT3TXaISNSMAM4Iw7XaLRzrvShUKhjtehUIi2Nu9zws997nN8/vOf56KLLuKZZ57h9ttvT2tNvuuTLy3MJWQ6kxeR4am+vp5x48YBsGzZso7lZ5xxBsuXLwdg5cqV7Nu3LyXH813Ih0NGmS6IEpFh6vbbb+eSSy7h1FNPpby8vGP5bbfdxsqVK5k+fTo///nPGT16NMXFxUM+njnnhryTVKmqqnJr1qwZ8n4W3v9HxpTk8dA/zE5BVSIyXG3cuJGpU6dmu4yUaG5uJhwOE4lEWLVqFZ/61KdYu3Ztj/WSfc9m9rJzrirZfn3XJw+66lVEjjzbtm3j0ksvJR6Pk5uby4MPPpiS/fo25De/sz/bZYiI9Nvxxx/PX/7yl5Tv13d98uCFfO2BZuLx4dMVJSLpMZy6nNNtMN+rP0O+KEprzFHf2JrtUkQkjfLy8ti7d28ggt45x969e8nLyxvQdr7trgFvrPyowtwsVyMi6TJ+/Hiqq6upqanJdikZkZeXx/jx4we0jb9Dfn8zk48e+hAkERmecnJyOq4YleR82V1TrknKREQAn4a85q8REfH4MuRH5EXIjYQU8iISeL4MeTOjoijKHoW8iAScL0MeoLwol70HW7JdhohIVvk25EcV5rJPIS8iAefbkC8tyKVOIS8iAZeSkDezpWa2x8w2dFlWamZPmtmWxGNqbnPST6MKc9nXoJAXkWBL1Zn8D4AFhyy7GXjKOXc88FTidcaUFubS0BKjqTWWycOKiAwrKQl559yzQN0hiy8G2m97sgz4UCqO1V+jCrzpDHQ2LyJBls4++aOdc7sSz3cDRydbycyuNbM1ZrYmlfNPlBbmAKhfXkQCLSMfvDpvirik08Q555Y456qcc1UVFRUpO2bHmfxBzUQpIsGVzpB/x8zGACQe96TxWD2UFqq7RkQknSH/GHBV4vlVwK/SeKweRinkRURSNoTyJ8AqYIqZVZvZ1cBdwPlmtgU4L/E6Y0bmq09eRCQl88k7567opencVOx/MCLhECX5ObrqVUQCzbdXvILXL1/XoA9eRSS4fB3yowp0Ji8iwebrkC8t1Pw1IhJsvg75UQWav0ZEgs3XId9+Ju9diyUiEjy+DvmRBbk0t8Vp1CRlIhJQvg55zV8jIkHn65DX/DUiEnS+Dvn2+Wvq9OGriASUr0O+Y/4addeISED5OuRLE9016pMXkaDydciPyM8hZJqJUkSCy9chHw4ZIwtyqd7XmO1SRESywtchD5CfE+bRv+zgpa2H3oJWRMT/fB3ybbE4NQeaAbhx+TraYvEsVyQiklm+DvllL2ztCPbaA838cNXW7BYkIpJhvg35PfubuPfJzcQT09Y0tMS4Z+VmavY3Z7cwEZEM8m3IP7Z2J7F494nJYnHHY+t2ZqkiEZHM823IXzxrHOGQdVsWDhkXzRybpYpERDLPtyFfURzlxvMnkxP2gj4vEuKm+ZOpKI5muTIRkczxbcgDXDWvsmOSshH5OVx5emV2CxIRyTBfh3wkHOKLC04A4OOnTSQS9vW3KyLSg+9T7+wpFYB3AxERkaDxfciPKsglZN44eRGRoPF9yIdDRmlhLrUHNEmZiASP70MeoKwwyl6dyYtIAAUi5MuLc9VdIyKBFIiQLyuMslc3DhGRAApEyJcXRanVnDUiEkCBCPmyolwOtsRobIlluxQRkYwKRMiXF3lj5Pce1Nm8iARLJN0HMLOtwH4gBrQ556rSfcxDlRd589XUHmhh/KiCTB9eRCRr0h7yCR9wztVm6Fg9lCVCXsMoRSRoAtVdo2GUIhI0mQh5B6w0s5fN7NoMHK+Hrt01IiJBkonumjOdczvM7CjgSTPb5Jx7tr0xEfzXAkycODEtBeTlhCmKRtirkBeRgEn7mbxzbkficQ/wKDDnkPYlzrkq51xVRUVF2uooK9JVryISPGkNeTMrNLPi9ufAfGBDOo/Zm7LCXA2hFJHASXd3zdHAo2bWfqwfO+eeSPMxkyovivL23oZsHFpEJGvSGvLOuTeBmek8Rn+VFUX587Z92S5DRCSjAjGEEqCiKJe6gy3E4i7bpYiIZExgQr6sKErcwb4GjbARkeAITMh3jpXXh68iEhyBCfnRJXkA7KpvynIlIiKZE5iQHzcyH4Cd7zZmuRIRkcwJTMhXFEeJhEwhLyKBEpiQD4eM0SV57HxX3TUiEhyBCXmAsSX57NCZvIgESLBCfmSeumtEJFACFvL57K5v0gVRIhIYgQv5trijZr/GyotIMAQq5NuHUapfXkSCIlAhP1Zj5UUkYAIW8t5Vrwp5EQmKQIV8cV4OxXkRhbyIBEagQh68fvkduiBKRAIicCE/dmS+zuRFJDACGPJ57KxXyItIMAQw5PN5t6GVhpa2bJciIpJ2gQv5zimH1S8vIv4XuJDXWHkRCZLAhryuehWRIAhcyI8ekUc0EuKNPQeyXYqISNoFLuTDIeP4o4t47Z392S5FRCTtAhfyAFOOHsGm3Qp5EfG/QIb8CaOLqdnfTN3BlmyXIiKSVoEM+SmjiwHYtPu9LFciIpJegQz5ExIhv1ldNiLic4EM+YriKKMKcvThq4j4XiBD3syYMrpYH76KiO8FMuQBThg9gs279+vKVxkWdtf3Ps1Gqtv62kb8J7AhP2V0MQdbYsy76w+8tLUu2+VIF8Ml8DLVtvqtOk6/86mkP4epbutrG/GntIe8mS0ws9fM7HUzuznlB1i/HO6bDreP9B7XL+/XZu+rKOx4fuPydbTF4ikvzU+SBlfivXfJ3vvBtK1fTvPdJ3LU14+m+e4Te2yT0rZMHquPttjanzHxh3N4I/oxJiybQ2ztzzo2SXVbX9uIf5lzLn07NwsDm4HzgWrgJeAK59yrydavqqpya9as6f8B1i+HFddDa5cul5x8+OADcNKlfW767ae3sPh3mwEoyA1z0/zJ/OOZx/b/2Eei9cvhqTtw9dVYyXg499bO96m3tvXLaf7d7eQc2EFr0TiiF9zesbzX9x4G3jbzY7h1P8a6LHc5+VhiG7fi+tS1ZfJYh6mj9c+PkBPv/APaGsoj50Pf9J7//8+lru2Uj/d+rMP8rsjwZ2YvO+eqkralOeRPB253zl2QeH0LgHPuzmTrDzjk75sO9dt7Li+ZADds6HWzPfubOPvuZ2hoiXUsK8gN8z9f+AAVxdH+Hz9bEoFMfTUMIKwHFLyRPJhxCe6Vn2NtncHgIlHszBvgpYfgYE3P2vJLAQeN+3q2RUu8tuZk1yeY13aonALvsbWhZ1vEm2yOtiSfq0QS/45tzUmO1YtwrvcYa6H916K9ItfR1ta5rL3ucMRbP96W+D66fCehiPc63vmzdqjO7zrxHlg40RDvuY6FerTRo83R87306nrHlXC01XdZbJ3v8SHVdO6r76oPt67hkq2NJd1H8r1bL1UkrSPpsuTNSWvoh/5sNdB9v37cVUz/xOJB1dNXyEcGtcf+Gwd0TeFq4LSuK5jZtcC1ABMnThzY3uurB7Y84bG1O4nFu/8DxOKOx9bt5OozJw2shnRKFubQPZDrt3uvwfvFX7EI2pq8X4r67fCrT8MbT8Gm33YPcfBe/+rT3i9n/JCbqLQ1wV9+1OOXy9qa4Zm7eq+5sY++3ubOcHHO+0VxWJcvcISIJx4dEG9xxDEchYlHI06IOEa81XDOiFHYsTxGyHveYt46iXW99vbXIZxrf91lnzHvefs+OmsyiHd53u0R3CFtXX+9Xbzn+u2vD13GoftOmoHJ99d1Wefr7vsFqHEjKLfuo8rcIf/0PfdjWJdF1vWZtb+2Q5LYvJfWfauuy6zLsfra/yFbJ1nRkry0JKslWXbowXvffY89Jt3qkHVdsnV77M9TPGJG0jqGKt0hf1jOuSXAEvDO5Ae0ccn4Xs7kx/e52cWzxvH1Jzd3WxYOGRfNHDugw6fVoWfe9dvhsc9BKCd5WD96XfKzu1grrPtp78eJtfZY1OaMZnJpdLk0EqWRKA1EaXB5HCTKQfI5GB5FQ5vjgMvjAPk0kOc9txEcIMrBWA6NRNlPPo3kESNEjDCxRMgORg5tRGgjlzZyiBE2R8RiRF0zEWKJ9hgRYlgoTNji5MUbiBAnQhthHBHaCBEnYnFyEtVEEo/xSD4h4hS0vUeYGGGLE078WWiLFCXa6juWhRNfLTkjCFmcwtZ9WMdy789Uc24JIRzFLXs72gwSf8qMsMUSy1zizxFQUI4BOQ27Ov4UtK9D4VEARA7u6tzGEmfLhaMxHDkHd9L+Z9La/1xaiJCL0epCtBKhiEbC5ogXjyWy6OXEPizxCCFLPO8lAOXIke6Q3wFM6PJ6fGJZapx7a/IuiPYz3l5UFEe58fzJ3L3yNZpa40RCxk3zJ2enqyZZ98px58ITt/QM87YmoJeRG8kCHu9scD8F7MqfQt3BJmopYa8bwV5XQh3F1IWPos4VUt+Ww7uukHcpoonDvw85bXGK3QEKrYkCmiigmSgtNBVPIjcUZ2r9yxRZE3k0k28thF2M3WPPI2KOE3f+MtHWQi6tGI5d5fM4vvYpRlgDObQRtTZizrAFXyUnBPm/XURBqHOuoUaXS+vf3g9EyPnNjeRbsjbI+c2iHm2N0y4j/68/G9A2rQv7aLuwj7YLem9zsz5Ozis/6dlPvtAbn5C0b/3CL/fetuBfem/r6JM/2H35/H+FnDDiX+kO+ZeA481sEl64Xw58LGV7b++HfvJW2L8L8kfBhYv79UHSVfMqWbbqbbbVNRAy48rTK1NWVr91OVvv6F7572vpX4+fp9lF2OHK2V4wjW1tI9ndGGanK2OnK2M3pbzjRnmhfUjXdIg4JRwgPyePwnCcibFXmRF6k5EcoMQOkk8zB0aewOR3n6c89B5FNFBsjURcGzkL7yQ+7aMsvvsOFvFTxtpedroy7udyvvjZiwBYfPd6/olfdm+78rJE2+ae2/3jrWx4vJKKv97HGPayizK2zLiBs+d5Py5Pb6tjcte26Tdw9pxE29YBtl3yWZ5efkLq9jfYtg9/ltik09n9qy9xVLyWPaFyKi76SsfPbyjuUtc26zJC40/rdRvxr7R+8ApgZguBbwBhYKlz7iu9rTvgD17bNb0Hd02A+f8B8z7X781Wv1XHpd9b5T3/0rkcVZw38GP3x6Fn6+f8Gxw1FX6wEJqTXHWbV+J9CNjlg81aN4LN8fFsDk9mc0sp2zmKN+Nj2EVZt+6PCG0czT7G2l5GWx1l9h6NY+Yw7dS/4anf/IyreYypoW00uCj/aZfxxS94/+tJGthfuJUNjz+YNCQBHvrjm9zz5GYaW2IU5IS56YLOEUqDaWuLxTnn3v9hW10DE0sL+MON7ycS9r63VLdl8liHa1v9Vh2XfW8Vy//5dGZXlnb7UUh1W1/byJEra6NrBmrQIR+Pwx2lcNYX4JwvDWjTP26u4ZNLV3P335/EJVUTDr/BQCUb1dLbSBK87pVqV8GqWV/ltZef4TU3kU3xidRS0rFOhDbeZzuYbNWU8R6v2PuYNuc8/r5qAsu+u5h/sZ5hXVEcHfbBC8Mn8DLdtru+idElyU8yUt3W1zZyZPJ/yAPcORFmXQEXfm1AmznnmHvnU1QdU8q3Pn7K4I7dl96GeeaNIhbOgwO72eAmsTp+Ai/Fp/Dn+OSOQDfiHG87OMneYKzV8WemMHXeB3n4xa00tnb2wXcd/tlXkB8JwQvDJ/Ay3SYyWMEI+fumQ+XfwIe/M+BNb/7lelas28kLN59LSUHO4I6f0O2XuKke7uo5LHR7vIJn4rP4ZeGlvFFv7Mcbp1xpu5lprzNi/FTOP+9CrvvRmqRh/tjaHYMKcjgygldEBqavkPfP3DXREb1cZHN4V82rpKE1xpI/vjH443e9bH3xVPjNTfCfs73x0w5ejU9kcetlnN+8mL9puZ9/a/s/7KGM8tB73B7+AS/mfoaHc75KbiTE9Z+8jLMmV3DT/Cnk53ojHwpywh0jgK6aV0lFkTcCprw42u1D40g4xD2XzMSAey+d2S3gAeZMKmXVLecm7Y/tq62vME5Hm4ikhn9CPm+Ed+Y8CFPHjODvThrLfz2/ldoDncNQ+j2p1frluBXXEz24g5BBtGEnvPQgjZbHnYX/lzOaH2Bhy118L/Z3lFs9N4d/zJ3Tqnn+5nP4xIUf4Gvhq5nb8i3mu28xdf7VHUM5ewvzoQQ5KHhFgsRHIV8y6JAHuOG842lui/Ptp72z+QHN4vfUHR3zkrS5EHEH/9V2AdNqvsJD+2YRtzD/Gn6Yl3I/xd2R7/FGuJLzPvQPmNmgz8qHEuQiEhz+CfkhdNcAHFtRxEdPGcfDL77Nuu37uPHn63B0n6GyLRZPurzp3d28EJvKQReljmIub/k3fh2byx2RH7DmS+fxTwvncV/4Kk5pWdLjbH0oZ+UKchE5nKxPa5AyeSO88fJD8Pnzp/Dcllou/d6LHQMcaw8088NVW/nHM49l2Qud3Tk1+5v46m830hZ31LR+mq/nfIdaV8Ljsdl8LedBJoV2sz9vDMWFud0uvDr0bB06g7y30FaYi8hg+e9MfgijhUaX5PHdT55KS1ucljbvLL2hJcY9Kzfz6s567ln5Go2JmSsbW+MsfX4rm176A/fnfptWwoy1Gq7N+S2TQru9y9bP8S40OtzZevuxRURSzT8hnzfCm0kx2ZS0A7D6rTpywt0nZWpoibHwgee6DWcEOC70Dsvy7yV31HjWT7mB3ZQTd8YOV86fpt/GiDmdMzgcrg9dRCQdfBTyiStCh9hlc/GscT3OtCMh4+NzJnYL/yIaWJrzNaKRMHzil8y97AtcUfh9jm1+hCsKv8+ZH/l0j33rbF1EMs0/IR8d4T0O4cNX6Jyhsuv49H9deAJf+cgMlsx8k+ej1/Nm9GOsjn6GibaH0OWPQNlx/eqSERHJNP8kUYrO5KGX8enrl3P2lv9gnNUSMiiwZgiFuk1ZoC4ZERlu/BPyHWfygx8r3y7pWXmXsfDtLN7m3bmpC3XJiMhw4q8hlDCkC6K66jGscZC3GhQRySb/ncmnoLumXbez8t5uKXiYWw2KiGSTf0K+vU9+iB+89uqc/9fzhr/9uNWgiEg2+SfkcwvBwik9k+8mFPEutMovBQxKJsAHH9Dt00RkWPNPn7wZRIvTcyYfj8Oz98BR0+Cfn/NG1YiIHAH8lVZDmG64T68/CTUb4YxFCngROaL4K7GiJenprnn+ARgxHqZ/JPX7FhFJI3+FfF5J6rtrql+Gt5+D0z8N4aHdGlBEJNN8FvJDn264hxce8P6HcMqVqd2viEgG+CvkoyNScsVrh4N7YdNv4ORPeB/qiogcYfwV8qk+k9/wC4i3wqyPHX5dEZFhyF8hn4Ibh3Sz9scwegaMnp6a/YmIZJi/Qj5vBLg4tBwY2n7WL4d7psCutfDudu+1iMgRyD8XQ0H36YYH24e+fjmsuB7aZ5xsetd7Dbq6VUSOOP46k0/FjUOeuqMz4Nu1NvaYUlhE5Ejgr5DPS8FMlJpSWER8xF8hH23vrhnCMEpNKSwiPuKvkM9LQXfNubcCmlJYRPwhbSFvZreb2Q4zW5v4WpiuY3Xo+OD13cHvY9ypgIO8kWhKYRE50qV7dM19zrl70nyMTvmjvMeGfYPfx+YnvMfr/gdGVQ65JBGRbPJXd00k6o2waagd/D5eexyOOlEBLyK+kO6Q/6yZrTezpWY2KtkKZnatma0xszU1NTVDP2JBGRwcZMg31MHbL8CUC4deh4jIMDCkkDez35vZhiRfFwPfAY4DZgG7gHuT7cM5t8Q5V+Wcq6qoqBhKOZ7CCjg4yD8Wbz4NLgaTFfIi4g9D6pN3zp3Xn/XM7EHg10M5Vr8VlsO72wa37VvPet09Y09ObU0iIlmSztE1Y7q8/DCwIV3H6mYo3TVvPQvHnAFhf832ICLBlc40W2xmswAHbAWuS+OxOhVWeB+8Oufd3Lu/6quh7k2YfU36ahMRybC0hbxz7pPp2nefCssh3uaNlc9P+llvcm/90XucdFZayhIRyQZ/DaEEKCj3Hg/uHdh2W/8I+aVw1LTU1yQikiX+C/nC9pAfwAgb57z++MozIeS/t0REgst/idYe8gO5IGrfVqjfrq4aEfEd/4V8R3fNAEL+rWe9R4W8iPiM/0K+cBAh//bzUHgUlE9OT00iIlniv5AfzPw12/8EE08b2JBLEZEjgP9CHgZ2QdT+d7w++QmnpbUkEZFs8GfID2T+mu0veo8T5qavHhGRLPFpyJdDQz/HyW/7E0TyYMzM9NYkIpIF/gz5gXTXbH8Rxp4Ckdz01iQikgX+DPmu89f0pbURdq3zPnQVEfEhn4Z8l/lrerN+OXxjhrfen3/kvRYR8Rl/hvzh5q9ZvxxWXN/54WxDrfdaQS8iPuPPkC8s8x57G2Hz1B1eV01XrY3echERH/FpyCduI9jbBVH11QNbLiJyhPJnyB9u/pqS8QNbLiJyhPJnyB9u/ppzb4VwTvdlOfnechERH/FnyB9u/pqTLoVxVWCJb79kAnzwAW+5iIiP+PeO1QVlfU9tsH8XTFkIlz+SuZpERDLMn2fyAKXHQu3m5G37d3uTkk08PaMliYhkmn9DfvQM2LMJ2lp6tm1LTEo2UZOSiYi/+Tvk461Q+1rPtu1/gkg+jD4p83WJiGSQj0M+EeC7N/Rs27YKxp2qSclExPf8G/Jlx3ln67tf6b68eT/sWq9JyUQkEPwb8qEwHD0Ndq/vvnzLSnAxeN/52alLRCSD/Bvy4PXL736l+5TDG1d4N+2eMCd7dYmIZIjPQ366N91w+5w0rU2weSWcsNA70xcR8Tmfh3z7h6+Jfvk3n4bWgzD1g9mrSUQkg/x7xStAzSbv8adXeFMXjJwI0RKoPCu7dYmIZIh/Q379cnj8i52v67d7XxNO09BJEQkM/3bXJLsxCMDe1zNfi4hIlgwp5M3sEjP7q5nFzazqkLZbzOx1M3vNzC4YWpmD0NsNQBp6uSWgiIgPDfVMfgPwEeDZrgvN7ETgcmAasAD4tplldjhLrzcGmZDRMkREsmlIIe+c2+icSzI5DBcDP3XONTvn3gJeBzI7MP3cW70bgXSlG4OISMCkq09+HLC9y+vqxLIezOxaM1tjZmtqavqY/32gTrrUuxFIyQTAdGMQEQmkw46uMbPfA6OTNH3JOferoRbgnFsCLAGoqqpyh1l9YE66VKEuIoF22JB3zp03iP3uALp2fo9PLBMRkQxKV3fNY8DlZhY1s0nA8cDqNB1LRER6MdQhlB82s2rgdOA3ZvY7AOfcX4HlwKvAE8BnnHOxoRYrIiIDM6QrXp1zjwKP9tL2FeArQ9m/iIgMjX+veBUREcy51A5oGQozqwHeHuTm5UBtCss50un96E7vRye9F9354f04xjlXkaxhWIX8UJjZGudc1eHXDAa9H93p/eik96I7v78f6q4REfExhbyIiI/5KeSXZLuAYUbvR3d6PzrpvejO1++Hb/rkRUSkJz+dyYuIyCEU8iIiPuaLkDezBYk7UL1uZjdnu55sMrOlZrbHzDZku5ZsM7MJZva0mb2auIPZomzXlE1mlmdmq81sXeL9+HK2a8o2Mwub2V/M7NfZriVdjviQT9xx6lvAhcCJwBWJO1MF1Q/w7sYl0Abc6Jw7EZgLfCbgPxvNwDnOuZnALGCBmc3NbklZtwjYmO0i0umID3m8O0697px70znXAvwU785UgeScexaoy3Ydw4Fzbpdz7s+J5/vxfpmT3rwmCJznQOJlTuIrsCMvzGw88LfA97NdSzr5IeT7fRcqCS4zqwROBv6U5VKyKtE9sRbYAzzpnAvy+/EN4ItAPMt1pJUfQl6kT2ZWBPwS+Bfn3HvZriebnHMx59wsvBv5zDGz6VkuKSvM7O+APc65l7NdS7r5IeR1FyrplZnl4AX8I865/852PcOFc+5d4GmC+/nNGcBFZrYVr4v3HDN7OLslpYcfQv4l4Hgzm2RmucDleHemkoAzMwMeAjY6576e7XqyzcwqzGxk4nk+cD6wKatFZYlz7hbn3HjnXCVeZvzBOfeJLJeVFkd8yDvn2oDPAr/D+2BteeLOVIFkZj8BVgFTzKzazK7Odk1ZdAbwSbyztLWJr4XZLiqLxgBPm9l6vJOjJ51zvh06KB5NayAi4mNH/Jm8iIj0TiEvIuJjCnkRER9TyIuI+JhCXkTExxTyIiI+ppAXEfGx/wVeLNGTDntswgAAAABJRU5ErkJggg==\n", "text/plain": [ "
" ] }, "metadata": { "needs_background": "light" }, "output_type": "display_data" } ], "source": [ "fig, ax = plt.subplots(1,1)\n", "oqupy.helpers.plot_correlations_with_parameters(bath.correlations, tempo_parameters, ax=ax)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We could feed this object into the `oqupy.tempo_compute()` function to get the dynamics of the system. However, instead of that, we can split up the work that `oqupy.tempo_compute()` does into several steps, which allows us to resume a computation to get later system dynamics without having to start over. For this we start with creating a `Tempo` object:" ] }, { "cell_type": "code", "execution_count": 15, "metadata": {}, "outputs": [], "source": [ "tempo = oqupy.Tempo(system=system,\n", " bath=bath,\n", " parameters=tempo_parameters,\n", " initial_state=up_density_matrix,\n", " start_time=0.0)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We can start by computing the dynamics up to time $5.0\\,\\Omega^{-1}$," ] }, { "cell_type": "code", "execution_count": 16, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "--> TEMPO computation:\n", "100.0% 50 of 50 [########################################] 00:00:00\n", "Elapsed time: 0.6s\n" ] }, { "data": { "text/plain": [ "" ] }, "execution_count": 16, "metadata": {}, "output_type": "execute_result" } ], "source": [ "tempo.compute(end_time=5.0)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "then get and plot the dynamics of expecatation values," ] }, { "cell_type": "code", "execution_count": 17, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "" ] }, "execution_count": 17, "metadata": {}, "output_type": "execute_result" }, { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAkYAAAGwCAYAAABM/qr1AAAAOXRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjcuMSwgaHR0cHM6Ly9tYXRwbG90bGliLm9yZy/bCgiHAAAACXBIWXMAAA9hAAAPYQGoP6dpAABHo0lEQVR4nO3deVxVdeLG8efeC5d9ERFwQXHfBVfEtGXSTJ1S2xyn0qyscdrJqagZbZt00spKR80mtRzLpW3aLKNJa8INxH3NDZVVlE3Z7r2/PyCKnxsocLiXz/v1ugMczrn3Acfu4znf7/eYHA6HQwAAAJDZ6AAAAAD1BcUIAACgHMUIAACgHMUIAACgHMUIAACgHMUIAACgHMUIAACgnJvRAZyJ3W7X8ePH5efnJ5PJZHQcAABQBQ6HQ3l5eWrWrJnM5gufE6IYVcPx48cVHh5udAwAAHAJUlJS1KJFiwvuQzGqBj8/P0llv1h/f3+D0wAAgKrIzc1VeHh4xfv4hVCMquGXy2f+/v4UIwAAnExVhsEw+BoAAKAcxQgAAKAcxQgAAKAcY4wAAKgBNptNJSUlRsdosKxW60Wn4lcFxQgAgMvgcDiUlpamU6dOGR2lQTObzWrdurWsVutlPQ/FCACAy/BLKQoJCZG3tzcLABvglwWYU1NT1bJly8v6M6AYAQBwiWw2W0Upaty4sdFxGrQmTZro+PHjKi0tlbu7+yU/D4OvAQC4RL+MKfL29jY4CX65hGaz2S7reShGAABcJi6fGa+m/gwoRgAAAOWcvhjNmTNHERER8vT0VHR0tDZs2HDefRctWiSTyVTp4enpWYdpAQBAfebUxWjZsmWKjY3V1KlTlZSUpMjISA0dOlQZGRnnPcbf31+pqakVj8OHD9dhYgAAUJ85dTF69dVXNXHiRE2YMEFdunTRvHnz5O3trXfeeee8x5hMJoWFhVU8QkND6zDxuZXa7Np5PFeni0uNjgIAQK2ozhUeSZo7d6569OhRceP2mJgYffXVV7We02mn6xcXFysxMVFxcXEV28xmswYPHqyEhITzHpefn69WrVrJbrerV69eeumll9S1a9dz7ltUVKSioqKKr3Nzc2vuB/iNlJNnNPyNHyRJYf6eah3so9ZNfNQm2Kfs82AfhQd5y93i1D0WANBA/XKFZ968eYqOjtasWbM0dOhQ7dmzRyEhIec8pkWLFpo+fbrat28vh8OhxYsXa+TIkdq8efN537drgtMWo6ysLNlstrPO+ISGhmr37t3nPKZjx45655131KNHD+Xk5GjmzJkaMGCAduzYoRYtWpy1/7Rp0/Tcc8/VSv7fOpFfpEbe7jp5ukRpuYVKyy1UwoETlfaxmE1qGeSt/m2CNC4mQp2b+td6LgBA9TkcDp0pubwp45fKy91S7dlZGzZs0BNPPKH169erVatWWrJkiZKSkvT555/rP//5T43k+u0VHkmaN2+evvjiC73zzjt66qmnznnMDTfcUOnrv//975o7d67WrVtHMaopMTExiomJqfh6wIAB6ty5s+bPn68XXnjhrP3j4uIUGxtb8XVubq7Cw8NrPFefiCBtnnKdThYU6+CJAh3MLNDBrLLHgawCHcoq0JkSW8W29zekKLp1kO4aEKEhXULlxpkkAKg3zpTY1GXK14a89s7nh8rbWvW39nXr1umaa67R888/rwULFuiJJ57Q888/rx07dmjlypWV9n3ppZf00ksvXfj1d+5Uy5YtK2271Cs8v2Wz2bRixQoVFBRUeh+vDU5bjIKDg2WxWJSenl5pe3p6usLCwqr0HO7u7urZs6f2799/zu97eHjIw8PjsrNWVSMfqxr5WNWrZaNK2+12h9LzCrUnLU8rEo9q1fY0rT+YrfUHs9UswFN3xLTSH/q2VJDP5d0fBgDQsMTGxurWW2/VX/7yF0nS2LFjNXbsWI0cOVI9e/astO+f/vQn3XbbbRd8vmbNmp217VKu8Pxi27ZtiomJUWFhoXx9ffXxxx+rS5cuVfnRLpnTFiOr1arevXsrPj5eo0aNklR2r5T4+Hg9+OCDVXoOm82mbdu2afjw4bWY9PKZzSY1DfBS0wAvXd0xRKk5Z/TvdUe0dMMRHc8p1Mur9uj1b/dpZFQzjR8Qoa7NAoyODAANlpe7RTufH2rYa1fV0aNHlZCQoJkzZ1Zsc3Nzk8PhOOcwkqCgIAUFBdVIzqrq2LGjkpOTlZOTo5UrV2r8+PFas2ZNrZYjpy1GUlnTHT9+vPr06aN+/fpp1qxZKigoqLiGOW7cODVv3lzTpk2TJD3//PPq37+/2rVrp1OnTmnGjBk6fPiw7r33XiN/jGprGuClyUM76sHftdNnW45rccIhbT+Wq+Wbjmr5pqO6umMT/X10dzUP9DI6KgA0OCaTqVqXs4yya9cuSVKvXr0qtu3Zs0f9+vVT9+7dz9r/Ui+lXc4VHqvVqnbt2kmSevfurY0bN+r111/X/PnzL3jc5aj/f3IXMGbMGGVmZmrKlClKS0tTVFSUVq1aVXG67siRIzKbfx1/c/LkSU2cOFFpaWlq1KiRevfurZ9++qnWT8vVFk93i27tE65berdQ0pGTWvi/Q/pqe5q+35Opoa+t1TMjOusPfcNZqh4AcJacnBxZLL8O1s7OztbMmTMVGRl5zv0v9VJaTVzh+YXdbq80W7w2mBwOh6NWX8GF5ObmKiAgQDk5OfL3r5+zwn7OzNdfVmxR0pFTkqQrOzTR9Ju6qxlnjwCgxhUWFurgwYNq3bq1091JYf/+/Wrfvr2ef/553XrrrXrkkUeUlZWlY8eOVcxQqynLli3T+PHjNX/+/IorPMuXL9fu3bsrTmbMnj1bH3/8seLj4yWVTYAaNmyYWrZsqby8PC1dulT/+Mc/9PXXX2vIkCFnvcaF/iyq8/7NdCYX07aJr1b8aYCeGd5ZHm5mrd1bdvZo2cYjogMDAH7Rrl07Pf/883r99dfVs2dPNWvWTN98842aN2+u66+/vkZfa8yYMZo5c6amTJmiqKgoJScnV7rCI5UN0v75558rvs7IyNC4cePUsWNHXXvttdq4ceN5S1FN4oxRNTjDGaPf+jkzX5NXbNHm8rNHV3Voouk3d1fTAM4eAUBNcOYzRq6GM0a4qLZNfLXyTwP09PBOsrqZtWZvpq57ba2Wb0rh7BEAAOdAMXJxFrNJ913ZVl8+PEhR4YHKKyzVEyu36r73ElVQxL3ZAAD4LYpRA9EuxFcfThqguGFlZ49W70zXrfMSlJZTaHQ0AADqDYpRA2Ixm3T/VW21/P4YBftatTM1V6P/+T/tSq2dm+MCQEPB8ATj1dSfAcWoAYoKD9THf75CbZv4KDWnULfOS9CavZlGxwIAp+Pu7i5JOn36tMFJUFxcLEmyWKq++ve5MCutGpxtVtrF5Jwu0f1LNmndgWxZzCa9MLKb/hjd8uIHAgAqpKam6tSpUwoJCZG3tzeL6hrAbrfr+PHjcnd3V8uWLc/6M6jO+zfFqBpcrRhJUnGpXU99uFUfbT4mSZp0dVv95bqOMpv5iw0AVeFwOJSWlqZTp04ZHaVBM5vNat26tazWs2+oXp33b6e+JQgun9XNrFdui1R4kLdej9+nud//rJTs05p5a6Q8q3EzQgBoqEwmk5o2baqQkBCVlJQYHafBslqtlW4DdqkoRpDJZNJjQzooPMhbcR9t1edbU5WWU6i3xvVRkM/ZzRsAcDaLxXLZ41tgPAZfo8ItvVto8d395Ofppk2HT+qWuT8pI5fp/ACAhoNihEoGtA3Wx38eoOaBXjqQVaA/vr1emXm1eydjAADqC4oRztIuxE/vT+yvpgGe2p+RrzveXq/sgmKjYwEAUOsoRjinlo29tXRif4X4eWhPep7ueHu9Tp2mHAEAXBvFCOfVOthH79/XX8G+HtqZmqs7/7VBOWeYcQEAcF0UI1xQ2ya+en9itBr7WLXtWI7Gv7NBeYWUIwCAa6IY4aLah/ppyb3RCvR2V3LKKd21cKMKikqNjgUAQI2jGKFKOjf115J7ouXv6abEwyc1YdFGnS6mHAEAXAvFCFXWrXmA3rsnWn4ebtpwMFv3Lt6kwhKb0bEAAKgxFCNUS2R4oBbf008+Vot++vmE7n8vUSU2u9GxAACoERQjVFuvlo206O5+8nK3aM3eTD390TZxL2IAgCugGOGS9I0I0j9v7yWzSVqReFRvfrff6EgAAFw2ihEu2TWdQvTCqG6SpFdX79XKxKMGJwIA4PJQjHBZbo9upUlXt5UkPfXhVv1vf5bBiQAAuHQUI1y2v1zXUTdGNlOp3aE/vZeoPWl5RkcCAOCSUIxw2cxmk2bc2kP9Wgcpr6hUdy3coLScQqNjAQBQbRQj1AgPN4veurO32jbxUWpOoSYs2qh8VscGADgZihFqTKC3VYsm9FOwr4d2pebqz/9OYo0jAIBToRihRoUHeeudu/rIy92itXsz9bdPtrPGEQDAaVCMUON6tAjUm2N7ymySPtiYojn/ZY0jAIBzoBihVgzuEqrnbuwqSZr5zV6t2p5mcCIAAC6OYoRac2dMhCZcESFJenx5svamM40fAFC/UYxQq54Z3lkD2jZWQbFNE9/dpJzTJUZHAgDgvChGqFVuFrNm/7GXWjTy0uETp/Xg+0my2RmMDQConyhGqHVBPla9dWfZTLUf9mXp5VW7jY4EAMA5UYxQJ7o089eMW3tIkuavPaBPk48ZnAgAgLNRjFBnft+jWcUNZ5/8cKu2H8sxOBEAAJVRjFCnJl/XUVd3bKLCErvufy9RWflFRkcCAKACxQh1ymI26fU/9FTrYB8dO3VGD3DbEABAPUIxQp0L8HLXgnG95evhpvUHs/Xi5zuNjgQAgCSKEQzSLsRPr42JkiQtTjis5RtTjA0EAIAoRjDQkC6hih3SQZL010+2a9tRBmMDAIxFMYKhHrymnYZ0CVWxza5J/05kZWwAgKEoRjCU2WzSzFsjFR7kpaMnz+jxFcmyszI2AMAgTl+M5syZo4iICHl6eio6OlobNmyo0nEffPCBTCaTRo0aVbsBcVEBXu6ae3tvWd3M+nZXht764YDRkQAADZRTF6Nly5YpNjZWU6dOVVJSkiIjIzV06FBlZGRc8LhDhw5p8uTJGjRoUB0lxcV0ax6gZ2/oKkma8fUerT9wwuBEAICGyKmL0auvvqqJEydqwoQJ6tKli+bNmydvb2+988475z3GZrPp9ttv13PPPac2bdpc8PmLioqUm5tb6YHaM7ZfuEb3bC6b3aEH39+sjLxCoyMBABoYpy1GxcXFSkxM1ODBgyu2mc1mDR48WAkJCec97vnnn1dISIjuueeei77GtGnTFBAQUPEIDw+vkew4N5PJpL+P7qYOob7KzCvSI+8ny8Z4IwBAHXLaYpSVlSWbzabQ0NBK20NDQ5WWlnbOY3788Uf961//0oIFC6r0GnFxccrJyal4pKSw1k5t87a66Z+395a31aKEAyf02uq9RkcCADQgTluMqisvL0933nmnFixYoODg4Cod4+HhIX9//0oP1L52Ib6afnMPSdLs/+7Xd7vTDU4EAGgo3IwOcKmCg4NlsViUnl75TTM9PV1hYWFn7f/zzz/r0KFDuuGGGyq22e1l9+hyc3PTnj171LZt29oNjSq7MbKZNh3K1rsJh/XYsi364uGBatHI2+hYAAAX57RnjKxWq3r37q34+PiKbXa7XfHx8YqJiTlr/06dOmnbtm1KTk6ueNx444265pprlJyczPiheuiZEZ0V2SJAOWdK9MDSzSoqtRkdCQDg4py2GElSbGysFixYoMWLF2vXrl2aNGmSCgoKNGHCBEnSuHHjFBcXJ0ny9PRUt27dKj0CAwPl5+enbt26yWq1Gvmj4Bw83Cya/cdeCvBy15aUU3rpi11GRwIAuDinvZQmSWPGjFFmZqamTJmitLQ0RUVFadWqVRUDso8cOSKz2am7X4MXHuSt18ZE6u5Fm7Q44bBi2jbW9d2aGh0LAOCiTA6Hg/nQVZSbm6uAgADl5OQwELuOTftql+avOSB/Tzd9+cggxhsBAKqsOu/fnE6BU5h8XUdFhQcqt7BUD7+/WSU2u9GRAAAuiGIEp+BuMevNsT3l5+GmpCOnNOtb1jcCANQ8ihGcRniQd8X6Rv/8/mf9uC/L4EQAAFdDMYJTGdGjqcb2aymHQ3psebIy84qMjgQAcCEUIzidKb/vUnE/tcdXbJGd+6kBAGoIxQhOx8tatr6Rp7tZa/dmasEPB4yOBABwERQjOKUOoX6aekNXSdKMr/do85GTBicCALgCihGc1h/6hmtEj6YqtTv08AeblVtYYnQkAICToxjBaZlMJk27qbtaNPJSSvYZxX20TaxXCgC4HBQjODV/T3e9Oban3MwmfbE1Vcs2phgdCQDgxChGcHo9WzbS5KEdJUnPfrZD+zPyDE4EAHBWFCO4hPsGtdGg9sEqLLHrofeTVVRqMzoSAMAJUYzgEsxmk165NVJBPlbtSs3VjFV7jI4EAHBCFCO4jBB/T/2j/JYhb/94UGv3ZhqcCADgbChGcClDuoTqjv4tJUmPr9ii7IJigxMBAJwJxQgu55nhXdQupOyWIU+s3MoUfgBAlVGM4HK8rBa98YeeslrM+nZXuv69/ojRkQAAToJiBJfUpZm/nri+bAr/i1/sZAo/AKBKKEZwWXdf0Zop/ACAaqEYwWUxhR8AUF0UI7i0EH9PzbiFKfwAgKqhGMHlXds5VHf2byWpbAr/ifwigxMBAOorihEahGdGdFb78in8T37IFH4AwLlRjNAgeLpb9MbYX6bwZ+j9DSlGRwIA1EMUIzQYnZv+OoX/hc936kBmvsGJAAD1DcUIDcrdV7TWFe0a60yJTY8tS1aJzW50JABAPUIxQoNiNps089ZI+Xu6acvRHL0Zv8/oSACAeoRihAanaYCXXrqpuyRp9n/3K/FwtsGJAAD1BcUIDdLvezTTTT2by+6QHl2WrPyiUqMjAQDqAYoRGqxnR3ZV80AvpWSf0XP/2WF0HABAPUAxQoPl7+mu18ZEyWSSViQe1VfbUo2OBAAwGMUIDVq/1kGadFVbSVLcx9uUnltocCIAgJEoRmjwHh3cQd2a++vU6RJNXrFFdjurYgNAQ0UxQoNndTNr1pie8nQ364d9WVqccMjoSAAAg1CMAEntQnz1zPDOkqRpX+3W3vQ8gxMBAIxAMQLK3dG/la7u2ETFpXY98kGyikptRkcCANQxihFQzmQy6eVbeijIx6pdqbl6dfVeoyMBAOoYxQj4jRA/T00rXxX7rbUHtO7ACYMTAQDqEsUI+H+Gdg3TbX1ayOGQHl++RbmFJUZHAgDUEYoRcA5Tbuiq8CAvHTt1Rs+yKjYANBgUI+AcfD3c9NptUTKbpI+SjulLVsUGgAaBYgScR5+IIE26umxV7KdZFRsAGgSKEXABj1z766rYf1m5VQ4Hq2IDgCujGAEXYHUz67XbouThZtbavZlasu6w0ZEAALXI6YvRnDlzFBERIU9PT0VHR2vDhg3n3fejjz5Snz59FBgYKB8fH0VFRem9996rw7RwRu1D/fTUsE6SpL9/uUs/Z+YbnAgAUFucuhgtW7ZMsbGxmjp1qpKSkhQZGamhQ4cqIyPjnPsHBQXpmWeeUUJCgrZu3aoJEyZowoQJ+vrrr+s4OZzN+JgIDWwXrMISux5blqwSm93oSACAWmByOPGgiejoaPXt21ezZ8+WJNntdoWHh+uhhx7SU089VaXn6NWrl0aMGKEXXnjhrO8VFRWpqKio4uvc3FyFh4crJydH/v7+NfNDwGmk5RRq6Ky1yjlTood/106x13U0OhIAoApyc3MVEBBQpfdvpz1jVFxcrMTERA0ePLhim9ls1uDBg5WQkHDR4x0Oh+Lj47Vnzx5deeWV59xn2rRpCggIqHiEh4fXWH44n7AAT704qpskafZ/9yvpyEmDEwEAaprTFqOsrCzZbDaFhoZW2h4aGqq0tLTzHpeTkyNfX19ZrVaNGDFCb775poYMGXLOfePi4pSTk1PxSElJqdGfAc7nhshmGhXVTHaH9NiyZBUUlRodCQBQg9yMDlDX/Pz8lJycrPz8fMXHxys2NlZt2rTR1Vdffda+Hh4e8vDwqPuQqNeeG9lNGw5m6/CJ03rxi10V91YDADg/pz1jFBwcLIvFovT09Erb09PTFRYWdt7jzGaz2rVrp6ioKD3++OO65ZZbNG3atNqOCxcS4OWumbdFSpLe33BE8bvSL3IEAMBZOG0xslqt6t27t+Lj4yu22e12xcfHKyYmpsrPY7fbKw2wBqpiQNtg3TuwtSTpyQ+3Kiuf/w8BgCtw2mIkSbGxsVqwYIEWL16sXbt2adKkSSooKNCECRMkSePGjVNcXFzF/tOmTdPq1at14MAB7dq1S6+88oree+893XHHHUb9CHBik4d2VMdQP2XlFyvuo22sig0ALsCpxxiNGTNGmZmZmjJlitLS0hQVFaVVq1ZVDMg+cuSIzOZfu19BQYH+/Oc/6+jRo/Ly8lKnTp20ZMkSjRkzxqgfAU7M092i18ZEaeScH7V6Z7pWbDqq2/oycxEAnJlTr2NU16qzDgIajnlrftb0r3bLx2rRV49cqZaNvY2OBAD4jQaxjhFQX0wc1Eb9IoJUUGzTY8uTZbPzbw0AcFYUI+AyWcwmvXJbpHw93JR4+KTmrfnZ6EgAgEtEMQJqQHiQt569sask6bXVe7X9WI7BiQAAl4JiBNSQm3s11/Vdw1Rqd+jRZckqLLEZHQkAUE0UI6CGmEwmvXRTdzXx89D+jHz9Y9VuoyMBAKqJYgTUoCAfq16+pYckaeH/DunHfVkGJwIAVAfFCKhh13QM0R39W0qSJq/YolOniw1OBACoKooRUAueHt5ZbYJ9lJZbqL9+sp1VsQHASVCMgFrgbXXTq2OiZDGb9PnWVH2afNzoSACAKqAYAbUkKjxQj1zbXpL0t0+26+jJ0wYnAgBcDMUIqEV/vrqterYMVF5RqR5fvoVVsQGgnqMYAbXIzWLWrDFR8rZatP5gtt7+4YDRkQAAF0AxAmpZq8Y+mvL7LpKkmd/s0c7juQYnAgCcD8UIqANj+oZrSJdQldgcenTZZlbFBoB6imIE1AGTyaTpN3VXsK+H9qbn6+VVe4yOBAA4B4oRUEca+3ro5Vu6S5Le+d9BVsUGgHqIYgTUod91CmVVbACoxyhGQB17ZniXilWxn2FVbACoVyhGQB3zslr02pgouZlN+mJrqj5JPmZ0JABAOYoRYIDI8EA9XL4q9pRPdiglm1WxAaA+oBgBBvnz1W3Vi1WxAaBeoRgBBilbFbunfD3ctOFQtuat+dnoSADQ4FGMAAO1bOytZ2/sKkl6bfVebUk5ZWwgAGjgKEaAwW7u1VwjejRVqd2hR5clq6Co1OhIANBgUYwAg5lMJr00qruaBnjqYFaBXvxip9GRAKDBohgB9UCAt7teuS1SJpP0/oYUfb0jzehIANAgUYyAemJA22DdN6iNJOmpD7cqPbfQ4EQA0PBQjIB6JPa6DurazF8nT5do8ootsjOFHwDqFMUIqEc83Cx6/Q9R8nAz64d9WVr00yGjIwFAg0IxAuqZdiF++uuIzpKk6at2a3darsGJAKDhoBgB9dAd/Vvpd51CVFxq1yPvJ6uwxGZ0JABoEChGQD1kMpn08i09FOxr1Z70PL28ao/RkQCgQaAYAfVUsK+HZtwSKUl6538HtWZvpsGJAMD1UYyAeuyaTiEaF9NKkvT48i3Kyi8yOBEAuDaKEVDPPT28szqG+ikrv0iTV2yRw8EUfgCoLRQjoJ7zdLfojbE9ZXUz6/s9mUzhB4BaRDECnEDHsF+n8E/7crd2pTKFHwBqA8UIcBJ39m+lazuFqNhm18Pvb9aZYqbwA0BNoxgBTuKXKfwhfh7al5GvF7/YaXQkAHA5FCPAiTT29dCrt0VJkv69/oi+3pFmbCAAcDGXXYyys7Nlt9trIguAKhjYPlj3X9lGkvTkh1uVllNocCIAcB2XVIx27typ6dOna8CAAWrSpIlCQkI0btw4ffjhhyooKKjpjAD+n8ev66huzf116nSJHluWLJudKfwAUBOqXIz27Nmjxx9/XO3bt1f//v21ceNG/elPf1J6erq+/PJLtWrVSs8//7yCg4M1bNgwzZ07tzZzAw2a1c2sN/7QU17uFiUcOKG31h4wOhIAuASTo4qrxS1cuFDr16/XyJEjde2118pqtZ5zv0OHDunTTz/VZ599pm+//bZGwxotNzdXAQEBysnJkb+/v9FxAC3fmKInPtwqN7NJKycNUFR4oNGRAKDeqc77d5XPGE2YMEHz5s3TsGHDzipFSUlJFZ9HRETokUceqbNSNGfOHEVERMjT01PR0dHasGHDefddsGCBBg0apEaNGqlRo0YaPHjwBfcH6rtb+7TQiB5NVWp36JEPNiu/qNToSADg1GpkVlq/fv0UGxtbaduXX35ZE099QcuWLVNsbKymTp2qpKQkRUZGaujQocrIyDjn/t9//73Gjh2r//73v0pISFB4eLiuu+46HTt2rNazArXBZDLppdHd1TzQS4dPnNZfP97GLUMA4DJU+VLahfTs2VMjR47U4cOHtXDhQklSr169Kp1Jqg3R0dHq27evZs+eLUmy2+0KDw/XQw89pKeeeuqix9tsNjVq1EizZ8/WuHHjLro/l9JQX206lK0xb62Tze7QjFt66NY+4UZHAoB6o1YupV2IyWTSs88+q8jISN1yyy0qKSmp9X+1FhcXKzExUYMHD67YZjabNXjwYCUkJFTpOU6fPq2SkhIFBQWd8/tFRUXKzc2t9ADqoz4RQYod0kGSNOXTHdqfkW9wIgBwTjVSjH5pX48++qhuuOEG3XjjjTpz5kxNPPV5ZWVlyWazKTQ0tNL20NBQpaVVbdG7J598Us2aNatUrn5r2rRpCggIqHiEh/OvcNRff7qqra5o11hnSmx6cGmSCku4ZQgAVNclFaP8/Mr/Gv3+++8rPh8/frzuu+++847zqS+mT5+uDz74QB9//LE8PT3PuU9cXJxycnIqHikpKXWcEqg6i9mk126LUmMfq3an5emlL3cZHQkAnM4lFaOAgAB9+OGH5/3+6NGjlZ2dfcmhqiI4OFgWi0Xp6emVtqenpyssLOyCx86cOVPTp0/XN998ox49epx3Pw8PD/n7+1d6APVZiL+nXrktUpL0bsJhrdrOLUMAoDouqRg5HA7Nnz9fV1xxhQYOHKhHH31UGzdurOlsF2S1WtW7d2/Fx8dXbLPb7YqPj1dMTMx5j3v55Zf1wgsvaNWqVerTp09dRAXq1NUdQ3T/VWW3DHli5RYdPXna4EQA4DwueYzR5s2b1atXLw0cOFA7duzQoEGDNHny5JrMdlGxsbFasGCBFi9erF27dmnSpEkqKCjQhAkTJEnjxo1TXFxcxf7/+Mc/9Le//U3vvPOOIiIilJaWprS0tLMuDQLObvJ1HRUVHqjcwlI98kGySmzczxAAqsLtUg9cunSphgwZUvH11q1bNXLkSDVv3lyPPfZYjYS7mDFjxigzM1NTpkxRWlqaoqKitGrVqooB2UeOHJHZ/Gv3mzt3roqLi3XLLbdUep6pU6fq2WefrZPMQF1wt5j15tieGv76D0o8fFKzvt2rvwztZHQsAKj3Lmkdo+DgYP3444/q1Knyf2i/+OILPfbYY9q7d2+NBaxPWMcIzuaLral6YGmSTCbpvbujNbB9sNGRAKDO1fo6RlFRURULOf5Wu3btdOTIkUt5SgC1YESPphrbr6UcDumx5cnKzCsyOhIA1GuXVIxefPFFvfHGG7rzzjuVkJCggoICZWRk6KWXXlLr1q1rOiOAyzDl913UIdRXmXlFenzFFtnt3DIEAM7nkopR//79tW7dOqWkpGjQoEHy9/dX06ZNtXLlSr3yyis1nRHAZfCyWjT7j73k6W7W2r2ZmrvmZ6MjAUC9ddn3SsvIyFBiYqLsdruio6MVHOy6YxgYYwRntmzjET354TaZTdL7E/sruk1joyMBQJ2olTFGaWlpKio6e3xCSEiIhg0bphEjRlQqRQcOHKhGZAC17bY+4bqpZ3PZHdJD729WVj7jjQDg/6tyMVq5cqWCgoI0evRoLVy4UJmZmWfts379ej399NPq2rWrIiMjazQogMtjMpn04uhuahfiq4y8Ij22LFk2xhsBQCXVupS2f/9+/ec//9Gnn36qdevWqW/fvho+fLgOHjyozz//XJI0YsQIjRw5UkOGDDnvPcicFZfS4Ar2pudp5Oz/6UyJTbFDOujha9sbHQkAalV13r8veYzRiRMn9Pnnn+vLL79URESERo4cqZiYGJlMpksK7QwoRnAVKxOPavKKLTKbpCX3RmtAW9cdGwgAdVKMGiKKEVzJX1Zs0YrEowr29dCXjwxUiJ9rneEFgF/U+gKPAJzf8yO7qWOon7Lyi/TI+4w3AgCJYgQ0WF5Wi+bc3kveVosSDpzQ6/H7jI4EAIajGAENWLsQX700ursk6c3v9umHfWfPNgWAhoRiBDRwo3o219h+4XI4pEc/SFZ6bqHRkQDAMBQjAJp6Q1d1buqvEwXFemjpZpXa7EZHAgBDUIwAyNPdojl/7Ckfq0UbDmVrxtd7jI4EAIagGAGQJLVp4qsZt5atWD9/7QGt2p5qcCIAqHsUIwAVhndvqnsHtpYkTV6xVQcy8w1OBAB1i2IEoJInh3VSv4gg5ReVatKSJJ0uLjU6EgDUGYoRgErcLWbN/mNPNfHz0J70PMV9tE0skA+goaAYAThLiL+n5vyxlyxmkz5NPq53Ew4bHQkA6gTFCMA59WsdpLhhnSRJL36xU4mHTxqcCABqH8UIwHndM7C1hncPU4nNoQf+naSs/CKjIwFAraIYATgvk8mkl2+JVNsmPkrLLdTD77P4IwDXRjECcEG+Hm6af2dveVst+unnE3pl9V6jIwFAraEYAbiodiF+evmWHpKkud//rG92pBmcCABqB8UIQJX8vkcz3X1F2eKPjy/fwuKPAFwSxQhAlcUN76S+EY2UV1Sq+95LVF5hidGRAKBGUYwAVJm7xaw5t/dSqL+H9mfkK3b5FtntLP4IwHVQjABUS4ifp+bf2UdWi1mrd6brje/2GR0JAGoMxQhAtUWFB+rF0d0kSbO+3afVO9MNTgQANYNiBOCS3NYnXHcNiJAkPbYsWfsz8owNBAA1gGIE4JI9M6KzolsHKb+oVBPfTVTOGQZjA3BuFCMAl8zdYtY/b++l5oFeOphVoEc/2Cwbg7EBODGKEYDL0tjXQ/Pv7C0PN7P+uydTr7EyNgAnRjECcNm6NQ/QP24uWxl79n/366ttqQYnAoBLQzECUCNG9WyueweWr4y9Yov2pDEYG4DzoRgBqDFPDeukK9o11ulimya+u0mnThcbHQkAqoViBKDGuFnMmj22l8KDvHQk+7QmLUlSic1udCwAqDKKEYAa1cjHqrfH9ZWP1aKEAyc05dMdcjiYqQbAOVCMANS4jmF+emNsT5lM0vsbjmjRT4eMjgQAVUIxAlArru0cqqeHdZYkvfD5Tn2/J8PgRABwcRQjALXm3kGtdVufFrI7pIeWbta+dGaqAajfKEYAao3JZNKLo7qrX+sg5RWV6p7Fm5RdwEw1APUXxQhArbK6mTXvjt6/mamWqOJSZqoBqJ+cuhjNmTNHERER8vT0VHR0tDZs2HDefXfs2KGbb75ZERERMplMmjVrVt0FBRq4IB+r/jW+r3w93LT+YLb+9sl2ZqoBqJecthgtW7ZMsbGxmjp1qpKSkhQZGamhQ4cqI+PcAzxPnz6tNm3aaPr06QoLC6vjtAA6hPrpzT/2lNkkLduUon/9eNDoSABwFqctRq+++qomTpyoCRMmqEuXLpo3b568vb31zjvvnHP/vn37asaMGfrDH/4gDw+POk4LQJKu6RiiZ0Z0kSS99OUufbc73eBEAFCZUxaj4uJiJSYmavDgwRXbzGazBg8erISEhBp7naKiIuXm5lZ6ALg8d18RobH9wmV3SA+/n6xdqfy9AlB/OGUxysrKks1mU2hoaKXtoaGhSktLq7HXmTZtmgICAioe4eHhNfbcQENlMpn03I3d1L9NkPKLSjVh4Ual5RQaHQsAJDlpMaorcXFxysnJqXikpKQYHQlwCVY3s+bf0Udtm/goLbdQExZtVF5hidGxAMA5i1FwcLAsFovS0yuPT0hPT6/RgdUeHh7y9/ev9ABQMwK83bVoQj8F+3poV2quHli6mRvOAjCcUxYjq9Wq3r17Kz4+vmKb3W5XfHy8YmJiDEwGoDrCg7z1zl195OVu0dq9mUzjB2A4pyxGkhQbG6sFCxZo8eLF2rVrlyZNmqSCggJNmDBBkjRu3DjFxcVV7F9cXKzk5GQlJyeruLhYx44dU3Jysvbv32/UjwBAUo8WgXpzbNk0/g82pmjOf/k7CcA4bkYHuFRjxoxRZmampkyZorS0NEVFRWnVqlUVA7KPHDkis/nX3nf8+HH17Nmz4uuZM2dq5syZuuqqq/T999/XdXwAvzG4S6ieu7Gr/vbpDs38Zq9aNPLWqJ7NjY4FoAEyOThvXWW5ubkKCAhQTk4O442AWvDSl7v01toDcreY9O7d0Ypp29joSABcQHXev532UhoA1/PU9Z00ontTldgcuv+9TdqXnmd0JAANDMUIQL1hNpv0ym2R6tOqkXILS3XXwo3KyGONIwB1h2IEoF7xdLdowbg+ah3so2OnzuieRZtUUFRqdCwADQTFCEC908jHqkUT+qqxj1XbjuXoT0sSVVzKGkcAah/FCEC91Kqxj/51V195Wy36YV+WHl+xRXY7c0UA1C6KEYB6Kyo8UPPv7C13i0mfbTmu5z7bwQKQAGoVxQhAvTaofRO9eluUTCZpccJhvfkdC0ACqD0UIwD13g2RzfTsDV0lSa+u3qsl6w4bnAiAq6IYAXAK4wdE6OFr20uS/vbpdn2xNdXgRABcEcUIgNN4bHB73R7dUg6H9OiyzfpxX5bRkQC4GIoRAKdhMpn0/MhuGt49rGJ17K1HTxkdC4ALoRgBcCoWs0mvjYnSFe0aq6DYprsWbtSBzHyjYwFwERQjAE7Hw82i+Xf2UffmAcouKNad/9qg1JwzRscC4AIoRgCckq+HmxZN6Ks25bcOuX3Beu6rBuCyUYwAOK3Gvh56795oNQ/00oGsAt3x9nplFxQbHQuAE6MYAXBqzQO9tHRitEL9PbQ3PV93vL1eOadLjI4FwElRjAA4vVaNffTve/sr2Neqnam5Gr9wg/KLSo2OBcAJUYwAuIR2Ib5acm+0Ar3dlZxySncv3KjTxZQjANVDMQLgMjqF+eu9u6Pl5+GmDYeydd+7iSossRkdC4AToRgBcCndWwRo0d195W216Mf9Wfrzv5NUXGo3OhYAJ0ExAuByercK0r/G95WHm1nf7c7QIx9sVqmNcgTg4ihGAFxSTNvGemtcH1ktZn21PU2Pr9gim91hdCwA9RzFCIDLuqpDE825vZfczCZ9mnxck1ds4cwRgAuiGAFwaUO6hOqNsT3lZjbp483H9NhyyhGA86MYAXB5w7s31ew/9pK7xaTPthzXQ+9vVgnlCMA5UIwANAjXdwvT3Nt7V4w5eoDZagDOgWIEoMEY3CVU88f1ltXNrG92pmvSkkQVlbLOEYBfUYwANCjXdAzRv8b3kYebWfG7M1gEEkAlFCMADc6g9k20cEJfeblbtGZvpu5dvElniilHAChGABqoAW2DtWjCrytk372Ie6sBoBgBaMCi2zTWe/f0k6+HmxIOnNBd72xUXmGJ0bEAGIhiBKBB690qSO/d009+nmU3nv3jgvU6kV9kdCwABqEYAWjwerZspPcn9leQj1XbjuXo1vkJOnbqjNGxABiAYgQAkro1D9CKP8WoeaCXDmQW6Ja5P2l/Rp7RsQDUMYoRAJRr28RXKyfFqF2Ir1JzCnXrvARtSTlldCwAdYhiBAC/0TTASyvuj1FkeKBOni7R2AXr9OO+LKNjAagjFCMA+H8a+Vi19N5oDWwXrNPFNt29aKO+2pZqdCwAdYBiBADn4OPhpn/d1UfDu4ep2GbXA0uT9P6GI0bHAlDLKEYAcB4ebha9ObaXxvZrKbtDivtom/75/X45HA6jowGoJRQjALgAi9mkl0Z305+vbitJennVHk39zw6V2uwGJwNQGyhGAHARJpNJT1zfSX8d0VmS9G7CYd3/XqIKiriFCOBqKEYAUEX3Dmqjubf3koebWfG7M3Tb/ASl5xYaHQtADaIYAUA1DOveVO/f11+NfazacTxXo+b8T7tSc42OBaCGUIwAoJp6tWykj/98hdo28alYCHLN3kyjYwGoAU5fjObMmaOIiAh5enoqOjpaGzZsuOD+K1asUKdOneTp6anu3bvryy+/rKOkAFxJy8be+mjSFerfJkj5RaW6e9FGLV3PdH7A2Tl1MVq2bJliY2M1depUJSUlKTIyUkOHDlVGRsY59//pp580duxY3XPPPdq8ebNGjRqlUaNGafv27XWcHIArCPB217t3R+umns1lszv09MfbNO2rXbLbmc4POCuTw4kX5IiOjlbfvn01e/ZsSZLdbld4eLgeeughPfXUU2ftP2bMGBUUFOjzzz+v2Na/f39FRUVp3rx5F3293NxcBQQEKCcnR/7+/jX3gwBwag6HQ6/H79Osb/dJkkZ0b6pXbouUp7vF4GQApOq9fzvtGaPi4mIlJiZq8ODBFdvMZrMGDx6shISEcx6TkJBQaX9JGjp06Hn3LyoqUm5ubqUHAPx/JpNJjw7uoNfGRMrdYtIX21J1y7yfdOzUGaOjAagmpy1GWVlZstlsCg0NrbQ9NDRUaWlp5zwmLS2tWvtPmzZNAQEBFY/w8PCaCQ/AJY3u2UJL7olWkI9V24/l6oY3f1TCzyeMjgWgGpy2GNWFuLg45eTkVDxSUlKMjgSgnotu01j/efAKdW3mr+yCYt3xr/Va+L+D3EYEcBJOW4yCg4NlsViUnp5eaXt6errCwsLOeUxYWFi19vfw8JC/v3+lBwBcTItG3vpw0gCNLh+U/dxnOzV5xVYVltiMjgbgIpy2GFmtVvXu3Vvx8fEV2+x2u+Lj4xUTE3POY2JiYirtL0mrV68+7/4AcKk83S169bZI/XVEZ1nMJn2YdFS3zU/QccYdAfWa0xYjSYqNjdWCBQu0ePFi7dq1S5MmTVJBQYEmTJggSRo3bpzi4uIq9n/kkUe0atUqvfLKK9q9e7eeffZZbdq0SQ8++KBRPwIAF2YymXTvoDZ69+5+CvR219ajObpx9o9af4BxR0B95dTFaMyYMZo5c6amTJmiqKgoJScna9WqVRUDrI8cOaLU1NSK/QcMGKClS5fqrbfeUmRkpFauXKlPPvlE3bp1M+pHANAAXNEuWJ89OFCdm/orK79Yt7+9Xu8mHGLcEVAPOfU6RnWNdYwAXI7TxaV6YuVWfb617B9so6Ka6cXR3eXr4WZwMsC1NYh1jADA2Xhb3fTm2J6KG9ZJFrNJnyQf1w1v/qgdx3OMjgagHMUIAOqQyWTS/Ve11Qf39VfTAE8dzCrQ6H/+pPe4tAbUCxQjADBA34ggffnwIF3bKUTFpXb97dMdemBpknLOlBgdDWjQKEYAYJBGPla9Pb6P/jqis9zMJn25LU2/f/MHbUk5ZXQ0oMGiGAGAgX6Z0r9y0gC1aOSllOwzumXeT3r7hwNcWgMMQDECgHogKjxQXzw8SNd3DVOJzaEXv9ilie9u0on8IqOjAQ0KxQgA6okAL3fNvaOXXhjZVVaLWd/uytDQWWu1emf6xQ8GUCMoRgBQj5hMJt0ZE6GP/jxA7UN8lZVfrInvbtLkFVuUW8jAbKC2UYwAoB7q1jxAnz00UPdf2UYmk7Qy8aiuf22t/rc/y+hogEujGAFAPeXpblHc8M5afn+MWgZ563hOoW5/e72mfrpdZ4ptRscDXBLFCADqub4RQfrqkUG6PbqlJGlxwmENf+MHJR4+aXAywPVQjADACfh4uOnvo7tr8d39FOZftmL2rfN+0j9W7VZhCWePgJpCMQIAJ3JVhyb6+tErNbpnc9kd0tzvf9bw13/QTz8z9gioCRQjAHAyAd7uem1MlObd0UtN/Dx0IKtAf1ywXo8v36LsgmKj4wFOjWIEAE7q+m5N9W3sVbqzfyuZTNKHSUf1u1e+1/JNKayaDVwiihEAOLEAL3e9MKqbPpw0QJ3C/HTqdImeWLlVY95ap/0ZeUbHA5wOxQgAXECvlo302UMD9fTwTvJyt2jDwWwNe/0HvfrNHgZnA9VAMQIAF+FuMeu+K9tqdeyV+l2nEJXYHHrju/26ftZafbc7nctrQBVQjADAxbRo5K1/je+jubf3Uqi/hw6dOK27F23SuHc2aHdartHxgHrN5OCfEFWWm5urgIAA5eTkyN/f3+g4AHBReYUlmv3f/Vr44yEV2+wym6QxfVvq8es6KNjXw+h4QJ2ozvs3xagaKEYAnNWRE6c17atd+mp7miTJ18NND1zTThOuiJCnu8XgdEDtohjVEooRAGe3/sAJvfjFLm07liNJatHIS3HDOmt49zCZTCaD0wG1g2JUSyhGAFyB3e7QR5uPacbXu5WeWyRJ6tOqkZ4a1kl9IoIMTgfUPIpRLaEYAXAlp4tLNW/NAb219mcVltglld1y5PHrOqhHi0BjwwE1iGJUSyhGAFxRas4ZvRG/T8s3HZXNXvaWMKRLqGKHdFDnpvy3Ds6PYlRLKEYAXNmhrAK9Eb9PnyQfU3k/0ogeTfXY4PZqF+JnbDjgMlCMagnFCEBDsD8jT699u09fbE2VJJlN0qio5npkcHu1auxjcDqg+ihGtYRiBKAh2Xk8V699u1erd6ZLkixmk0ZGNtN9V7VRpzD+GwjnQTGqJRQjAA3RlpRTemX1Xq3dm1mx7eqOTXT/lW3Vv00Q0/xR71GMagnFCEBDtiXllN5ae0BfbU+tGIMU2SJA91/VVkO7hslipiChfqIY1RKKEQCUDdJ++8cDWrHpqIpKy6b5RzT21sQr2+jmXi1YSRv1DsWollCMAOBXWflFevenQ1qccFg5Z0okScG+Vo3t11Jj+7VUs0AvgxMCZShGtYRiBABnKygq1fJNKXr7h4M6duqMpLKZbL/rFKo7Y1ppULtgmbnMBgNRjGoJxQgAzq/EZtc3O9K1ZN1hJRw4UbG9VWNv3R7dUrf2DlcjH6uBCdFQUYxqCcUIAKpmf0aelqw7og+TjiqvsFSSZHUz6/c9muqO/q3UMzyQ2WyoMxSjWkIxAoDqOV1cqv8kH9d76w5rx/Hciu3tQ3x1c+8WGhXVXGEBngYmRENAMaolFCMAuDQOh0PJKae0ZN0Rfb71eMVsNrNJuqJdsG7p3ULXdQmTl5UZbah5FKNaQjECgMuXW1iiL7em6sOko9p46GTFdl8PNw3vHqabe7VQ34ggBmyjxlCMagnFCABq1uETBfoo6Zg+2nxUKdlnKraHB3lpePemGtatqSJbBDAeCZeFYlRLKEYAUDvsdoc2HsrWR0nH9MW2VOUXlVZ8r2mAp4Z2DdPQrmHq1zqIFbZRbRSjWkIxAoDad6bYpvjd6Vq1PU3/3Z2hgmJbxfca+1g1pEuohnYL04C2jeXhxpgkXBzFqJZQjACgbhWW2PS//Vn6anuavt2VrlOnSyq+5+fhpoHtg3V1xya6skMTNQ1gpW2cG8WollCMAMA4JTa7NhzM1qrtafp6R5oy8ooqfb9jqJ+u6thEV3Vooj4RjTibhAoUo1pCMQKA+sFud2jL0VNaszdTa/ZmKjnllH77bublbtGAto11VccmGtC2sdo28WUAdwNGMaolFCMAqJ9OFhTrx/1Z+n5PWVHKyq98NinIx6q+EY3Ur3VjRbcOUuem/gzibkBcvhhlZ2froYce0meffSaz2aybb75Zr7/+unx9fc97zFtvvaWlS5cqKSlJeXl5OnnypAIDA6v1uhQjAKj/7HaHdqXlas3eTP2wN0ubU06qsMReaR8/Dzf1jmikfq2D1C8iSN2aB8jTnUtvrsrli9GwYcOUmpqq+fPnq6SkRBMmTFDfvn21dOnS8x4za9YsFRYWSpLi4uIoRgDQQBSX2rXtWI42HMzWhoMntOnQSeX9ZjkASXIzm9Q+1E89mgeoe4sA9WgRoI5hfoxTchEuXYx27dqlLl26aOPGjerTp48kadWqVRo+fLiOHj2qZs2aXfD477//Xtdcc02VilFRUZGKin49HZubm6vw8HCKEQA4MZvdoV2pueVFKVsbD2XrREHxWfu5W0zqFOav7i0C1L15WVFqH+IrP093A1LjclSnGLnVUaYak5CQoMDAwIpSJEmDBw+W2WzW+vXrNXr06Bp7rWnTpum5556rsecDABjPYjapW/MAdWseoLsHtpbD4VBqTqG2Hs3RtmOnyj/m6NTpEm07Vvb5bzUN8FT70LKS1CHUV+1C/NQ+1Ff+FCaX4HTFKC0tTSEhIZW2ubm5KSgoSGlpaTX6WnFxcYqNja34+pczRgAA12EymdQs0EvNAr10fbcwSWU3vT168oy2Hs3R1mOntONYrvam5ykjr0ipOYVKzSnU2r2ZlZ4nzN9TLRt7q0UjL4U38lZ4UPnnQd4K8/dksPc5OBwO5ReVKiu/WJl5RcrKL5Knu1m/6xRqWKZ6U4yeeuop/eMf/7jgPrt27aqjNGU8PDzk4eFRp68JADCeyWRSeFBZuRnRo2nF9pzTJdqXkad9Gfnam56n/eUf03OLlJZbqLTcQm04ePbzuZnLyld4kJdC/DzVxM9Dwb5WBft6lH9e9rGRt9WpC1Spza7cwlKdPF2sU6dLdOqXj2fKPs/KL1JmXtnHss+LVFRaeWB8ZHggxUiSHn/8cd11110X3KdNmzYKCwtTRkZGpe2lpaXKzs5WWFhYLSYEADR0Ad7u6hMRpD4RQZW255wp0c+Z+UrJPq2jJ88oJfu0Uk6WfX7s5BmV2h06kn1aR7JPX/D5zSYpyMdDAV5u8vN0l7+Xu/w83eTvWf51+Uc/Tzd5uFlkdTPL3WKS1c0sq8Vc/rW54utf/DKa2CHHbz4vO2NTYnOouNSuolJb2UebXUUldhXb7BXbTxfZlF9UqvyiUhWUf/z1c5vyi0qUc7pEuYWlZ/9QVeBjtSjYz0NNfD3UqanfJT1HTak3xahJkyZq0qTJRfeLiYnRqVOnlJiYqN69e0uSvvvuO9ntdkVHR9d2TAAAzhLg5a5eLRupV8tGZ33PZncoLbdQR8tLU2Z+kbLyiso+5hcpK69YmflFOnm6WHaHKs6mODM/DzcFeLurkbdVgd7uCvS2KtDLXY3//1kyXw8F+1nlba03daT+FKOq6ty5s66//npNnDhR8+bNU0lJiR588EH94Q9/qJiRduzYMV177bV699131a9fP0llY5PS0tK0f/9+SdK2bdvk5+enli1bKigo6LyvBwDA5bCYTWoe6KXmgV660D/fS212ZReUlaTcM6XKKyw7A5NXWKK8Sh9LlVtYoqJSu0rKz+oU//bz33yUJJPKLs2ZTJJJqlgB3FT+P7+cafJwKz/TVH626ZczUlY3s3ysFvl4uMnX002+VrdfP/co/9zDogCvsgIU4OUu99+crXI2TleMJOnf//63HnzwQV177bUVCzy+8cYbFd8vKSnRnj17dPr0r6cs582bV2mG2ZVXXilJWrhw4UUv4QEAUNvcLGaF+HsqxN/T6CgNmtOtY2QkFngEAMD5VOf923nPdQEAANQwihEAAEA5ihEAAEA5ihEAAEA5ihEAAEA5ihEAAEA5ihEAAEA5ihEAAEA5ihEAAEA5ihEAAEA5ihEAAEA5ihEAAEA5ihEAAEA5N6MDOBOHwyGp7C69AADAOfzyvv3L+/iFUIyqIS8vT5IUHh5ucBIAAFBdeXl5CggIuOA+JkdV6hMkSXa7XcePH5efn59MJlONPndubq7Cw8OVkpIif3//Gn1u/Irfc93g91w3+D3XHX7XdaO2fs8Oh0N5eXlq1qyZzOYLjyLijFE1mM1mtWjRolZfw9/fn790dYDfc93g91w3+D3XHX7XdaM2fs8XO1P0CwZfAwAAlKMYAQAAlKMY1RMeHh6aOnWqPDw8jI7i0vg91w1+z3WD33Pd4XddN+rD75nB1wAAAOU4YwQAAFCOYgQAAFCOYgQAAFCOYgQAAFCOYlQPzJkzRxEREfL09FR0dLQ2bNhgdCSXs3btWt1www1q1qyZTCaTPvnkE6MjuaRp06apb9++8vPzU0hIiEaNGqU9e/YYHcvlzJ07Vz169KhYBC8mJkZfffWV0bFc3vTp02UymfToo48aHcWlPPvsszKZTJUenTp1MiwPxchgy5YtU2xsrKZOnaqkpCRFRkZq6NChysjIMDqaSykoKFBkZKTmzJljdBSXtmbNGj3wwANat26dVq9erZKSEl133XUqKCgwOppLadGihaZPn67ExERt2rRJv/vd7zRy5Ejt2LHD6Ggua+PGjZo/f7569OhhdBSX1LVrV6WmplY8fvzxR8OyMF3fYNHR0erbt69mz54tqex+bOHh4XrooYf01FNPGZzONZlMJn388ccaNWqU0VFcXmZmpkJCQrRmzRpdeeWVRsdxaUFBQZoxY4buueceo6O4nPz8fPXq1Uv//Oc/9eKLLyoqKkqzZs0yOpbLePbZZ/XJJ58oOTnZ6CiSOGNkqOLiYiUmJmrw4MEV28xmswYPHqyEhAQDkwE1IycnR1LZmzZqh81m0wcffKCCggLFxMQYHcclPfDAAxoxYkSl/1ajZu3bt0/NmjVTmzZtdPvtt+vIkSOGZeEmsgbKysqSzWZTaGhope2hoaHavXu3QamAmmG32/Xoo4/qiiuuULdu3YyO43K2bdummJgYFRYWytfXVx9//LG6dOlidCyX88EHHygpKUkbN240OorLio6O1qJFi9SxY0elpqbqueee06BBg7R9+3b5+fnVeR6KEYBa8cADD2j79u2GjhVwZR07dlRycrJycnK0cuVKjR8/XmvWrKEc1aCUlBQ98sgjWr16tTw9PY2O47KGDRtW8XmPHj0UHR2tVq1aafny5YZcGqYYGSg4OFgWi0Xp6emVtqenpyssLMygVMDle/DBB/X5559r7dq1atGihdFxXJLValW7du0kSb1799bGjRv1+uuva/78+QYncx2JiYnKyMhQr169KrbZbDatXbtWs2fPVlFRkSwWi4EJXVNgYKA6dOig/fv3G/L6jDEykNVqVe/evRUfH1+xzW63Kz4+nrECcEoOh0MPPvigPv74Y3333Xdq3bq10ZEaDLvdrqKiIqNjuJRrr71W27ZtU3JycsWjT58+uv3225WcnEwpqiX5+fn6+eef1bRpU0NenzNGBouNjdX48ePVp08f9evXT7NmzVJBQYEmTJhgdDSXkp+fX+lfHwcPHlRycrKCgoLUsmVLA5O5lgceeEBLly7Vp59+Kj8/P6WlpUmSAgIC5OXlZXA61xEXF6dhw4apZcuWysvL09KlS/X999/r66+/NjqaS/Hz8ztrfJyPj48aN27MuLkaNHnyZN1www1q1aqVjh8/rqlTp8pisWjs2LGG5KEYGWzMmDHKzMzUlClTlJaWpqioKK1ateqsAdm4PJs2bdI111xT8XVsbKwkafz48Vq0aJFBqVzP3LlzJUlXX311pe0LFy7UXXfdVfeBXFRGRobGjRun1NRUBQQEqEePHvr66681ZMgQo6MB1Xb06FGNHTtWJ06cUJMmTTRw4ECtW7dOTZo0MSQP6xgBAACUY4wRAABAOYoRAABAOYoRAABAOYoRAABAOYoRAABAOYoRAABAOYoRAABAOYoRAABAOYoRAABAOYoRgAblscce00033XTO723dulU33XSTGjduLE9PT3Xt2lUzZsxQaWlpHacEYBSKEYAGZcOGDerTp89Z29euXav+/fvLy8tLn376qbZs2aInn3xSr776qm666SbZ7XYD0gKoa9wrDUCDUFxcLB8fn0pnf6Kjo7Vu3TrZbDZ16NBBMTExWrJkSaXjdu/erR49emju3Lm655576jo2gDpGMQLQINjtdm3atEnR0dFKTk5WaGioPD09FRgYqISEBA0YMEDJycmKjIw869jRo0eroKBA33zzjQHJAdQlLqUBaBDMZrOOHz+uxo0bKzIyUmFhYQoMDJQkHTx4UJLUvn37cx7bvn17HT58uK6iAjAQxQhAg7F58+ZznhHy9/eXJGVnZ5/zuJMnT1bsA8C1UYwANBjnu1QWExMjd3d3ffbZZ2d9z2az6euvv9bAgQPrIiIAg1GMADQY27ZtU1RU1FnbGzdurIcfflgvvviijh8/Xul7r732mrKzs/XYY4/VUUoARqIYAWgw7Ha79uzZo+PHjysnJ6die35+vh5++GFFRETommuuUVJSkiRpxowZevrpp/Xmm2/KarXKZrMZFR1AHWFWGoAGY8mSJXryySd1/PhxTZ48WTNmzJAkPfvss3ruuecq9hs/frwWLVokk8lU6fiDBw8qIiKiLiMDqGMUIwAAgHJcSgMAAChHMQIAAChHMQIAAChHMQIAAChHMQIAAChHMQIAAChHMQIAAChHMQIAAChHMQIAAChHMQIAAChHMQIAACj3fx0SmdE68az1AAAAAElFTkSuQmCC", "text/plain": [ "
" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "dynamics_2 = tempo.get_dynamics()\n", "plt.plot(*dynamics_2.expectations(0.5*sigma_z, real=True), label=r'$\\alpha=0.3$')\n", "plt.xlabel(r'$t\\,\\Omega$')\n", "plt.ylabel(r'$\\langle\\sigma_z\\rangle$')\n", "plt.legend()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "then continue the computation to $15.0\\,\\Omega^{-1}$," ] }, { "cell_type": "code", "execution_count": 18, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "--> TEMPO computation:\n", "100.0% 100 of 100 [########################################] 00:00:02\n", "Elapsed time: 2.4s\n" ] }, { "data": { "text/plain": [ "" ] }, "execution_count": 18, "metadata": {}, "output_type": "execute_result" } ], "source": [ "tempo.compute(end_time=15.0)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "and then again get and plot the dynamics of expecatation values." ] }, { "cell_type": "code", "execution_count": 19, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "" ] }, "execution_count": 19, "metadata": {}, "output_type": "execute_result" }, { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAkYAAAGwCAYAAABM/qr1AAAAOXRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjcuMSwgaHR0cHM6Ly9tYXRwbG90bGliLm9yZy/bCgiHAAAACXBIWXMAAA9hAAAPYQGoP6dpAABJZElEQVR4nO3deVzUdeI/8NcczAzXDHILgqDgraAgiEdmUmZmWVtZ26ZLrftry9JoW7NddbMDzWrtcNXctdrKb26HHW5pRmqZeIF4izciOBxyDAwyM8x8fn8wjJGIiDPzmeP1fDzmgXzmMzOvkWRevT/vz/sjEQRBABERERFBKnYAIiIiIlfBYkRERERkxWJEREREZMViRERERGTFYkRERERkxWJEREREZMViRERERGQlFzuAO7FYLCgrK0NgYCAkEonYcYiIiKgTBEFAfX09oqKiIJV2PCbEYnQNysrKEBMTI3YMIiIi6oKSkhL06NGjw31YjK5BYGAggJa/WLVaLXIaIiIi6gydToeYmBjb53hHWIyuQevhM7VazWJERETkZjozDYaTr4mIiIisWIyIiIiIrFiMiIiIiKw4x4iIiMgOzGYzTCaT2DG8lkKhuOqp+J3BYkRERHQdBEGAVqtFbW2t2FG8mlQqRXx8PBQKxXU9D4sRERHRdWgtReHh4fDz8+MCwCJoXYD5/PnziI2Nva6fAYsRERFRF5nNZlspCgkJETuOVwsLC0NZWRmam5vh4+PT5efh5GsiIqIuap1T5OfnJ3ISaj2EZjabr+t5WIyIiIiuEw+fic9ePwMWIyIiIiIrty9Gy5YtQ1xcHFQqFdLT07Fr164r7vvee+9BIpG0ualUKiemJSIiIlfm1sVo7dq1yM7OxoIFC1BQUICkpCRMmDABFRUVV3yMWq3G+fPnbbfi4mInJiYiIiJX5tbF6PXXX8eMGTOQlZWFAQMGYMWKFfDz88Pq1auv+BiJRILIyEjbLSIiwomJr2z7ySoYmq9vwhgREZGrupYjPACwfPlyDBkyxHbh9oyMDHz77bcOz+m2xchoNCI/Px+ZmZm2bVKpFJmZmcjLy7vi4xoaGtCzZ0/ExMTgzjvvxKFDh664r8FggE6na3NzhJLqRvzuXzsxatFmvJl7HBcaDA55HSIiIjF05QhPjx49sGjRIuTn52PPnj246aabrvq5bQ9uW4yqqqpgNpsvG/GJiIiAVqtt9zF9+/bF6tWr8eWXX+LDDz+ExWLByJEjce7cuXb3z8nJgUajsd1iYmLs/j6AlmIUHqhCVYMBr286hpGLfsBbuccd8lpERORYgiCg0dgsyk0QhGvOu2vXLtx4443w9fVFv379sGfPHrzzzju444477PZ30pUjPJMnT8Ztt92GxMRE9OnTBy+99BICAgKwY8cOu+Vqj1ct8JiRkYGMjAzb9yNHjkT//v2xcuVKvPDCC5ftP3fuXGRnZ9u+1+l0DilHIxNC8dOccfjmwHn866fTOFBah9c2HUNUkC9+k9LD7q9HRESOc9FkxoD5G0V57cMLJ8BP0fmP9h07dmDcuHFYuHAhVq1ahb/85S9YuHAhDh06hE8//bTNvi+//DJefvnljl//8GHExsa22dZ6hGfu3Lm2bZ05wvNLZrMZn3zyCfR6fZvPcUdw22IUGhoKmUyG8vLyNtvLy8sRGRnZqefw8fHB0KFDceLEiXbvVyqVUCqV1521U1lkUtyZHI07kqLwj03H8OYPJ/DcugPo312NAVFqp2QgIiLvkp2djXvvvRfPPPMMAOCBBx7AAw88gDvvvBNDhw5ts++jjz6K++67r8Pni4qKumxbR0d4jh492uHzHThwABkZGWhqakJAQADWrVuHAQMGdOatdZnbFiOFQoGUlBTk5uZiypQpAFqulZKbm4uZM2d26jnMZjMOHDiA2267zYFJr41EIsHszD7YX1qHLUWVePTDfHw9czQ0fl1f3pyIiJzH10eGwwsniPbanXXu3Dnk5eXh1VdftW2Ty+UQBAHPP//8ZfsHBwcjODjYLjk7q2/fvigsLERdXR0+/fRTTJ8+HVu3bnVoOXLbYgS0NN3p06cjNTUVaWlpWLp0KfR6PbKysgAA06ZNQ3R0NHJycgAACxcuxIgRI5CQkIDa2losWbIExcXF+MMf/iDm27iMVCrB0qnJuP2tbThb3Yg5n+3HiodSxI5FRESdIJFIrulwlliOHDkCABg2bJhtW1FREdLS0jB48ODL9u/qobTrOcKjUCiQkJAAAEhJScHu3bvxxhtvYOXKlR0+7nq4/k+uA1OnTkVlZSXmz58PrVaL5ORkbNiwwTZcd/bsWUill+aX19TUYMaMGdBqtejWrRtSUlKwfft2hw/LdUWQnwIrfpeCO5f9jA2HtCg4W4Nhsd3EjkVERB6irq4OMpnMdimN6upqvPrqq0hKSmp3/64eSrPHEZ5WFosFBoNjz9yWCF2Zwu6ldDodNBoN6urqoFY7Z97PXz7dh//uOYexfcLw/sNpTnlNIiLqnKamJpw+fRrx8fFudyWFEydOIDExEQsXLsS9996LWbNmoaqqCqWlpdi5cyd69uxpt9dau3Ytpk+fjpUrV9qO8Pz3v//F0aNHbYMZb7/9NtatW4fc3FwALSdATZw4EbGxsaivr8eaNWuwePFibNy4ETfffPNlr9HRz+JaPr/d9nR9bzFzXCJkUgm2HqvE3rM1YschIiIPkZCQgIULF+KNN97A0KFDERUVhe+++w7R0dG49dZb7fpaU6dOxauvvor58+cjOTkZhYWFbY7wAC2TtE+ePGn7vqKiAtOmTUPfvn0xfvx47N69+4qlyJ44YnQNxBgxAoBnPtmHT/LP4ca+YXgvi6NGRESuwp1HjDwNR4y8yMybEiCTSrClqBKFJbVixyEiIvJYLEZuoGeIP+4eGg0AeOP7YyKnISIi8lwsRm7i8XEJkEiAzUWVKKluFDsOERGRR2IxchNxof4YnRAKAPgkv/1ruxERkTg4XVd89voZsBi5kXtTW67T9umeEpgt/EdIRCQ2H5+WqxI0NnIkX2xGoxEAIJN1fvXv9rj1Ao/e5pYBEdD4+qCsrgnbTlRhbJ8wsSMREXk1mUyGoKAgVFRUAAD8/PxsCyaS81gsFlRWVsLPzw9y+fVVGxYjN6LykWFKchTezyvGf/eUsBgREbmA1statJYjEodUKkVsbOx1F1MWIzdz3/AYvJ9XjE2HylGjN6Kbv0LsSEREXk0ikaB79+4IDw+HyWQSO47XUigUbS4D1lUsRm5mYJQGA6PUOFSmw7q9pXh4dLzYkYiICC2H1a53fguJj5Ov3dDU4S2TsP+7p0TkJERERJ6FxcgN3ZEUBblUgqPaepysbBA7DhERkcdgMXJDQX4KjLSuabThoFbkNERERJ6DxchNTRzUchbEtwfPi5yEiIjIc7AYualbBkRAKgEOlupw9gIXFiMiIrIHFiM3FRKgRHp8CABgwyGOGhEREdkDi5Ebu21wy+G0bw5wnhEREZE9sBi5sQkDIyGRAIUltSirvSh2HCIiIrfHYuTGwtUqpPbsBoBnpxEREdkDi5Gbu3VQdwAsRkRERPbAYuTmbrWetr+nuBo1eqPIaYiIiNwbi5Gbiw7yRb/IQFgE4MfjlWLHISIicmssRh5gXL9wAMAPRytETkJEROTeWIw8wLi+LcVo67FKmC2CyGmIiIjcF4uRBxgWGwSNrw9qG00oLKkROw4REZHbYjHyAHKZFDf0CQMAbD7KeUZERERdxWLkIcb1bSlGnGdERETUdSxGHmJsnzBIJMDh8zpo65rEjkNEROSWWIw8REiAEkk9ggAAW4o4akRERNQVLEYe5Caetk9ERHRdWIw8SGsx+vlEFUxmi8hpiIiI3A+LkQcZ0F2NEH8F9EYz9p6tFTsOERGR22Ex8iBSqQQjE0IBANt4eRAiIqJrxmLkYca0FqMTVSInISIicj8sRh5mVGJLMdp3rg66JpPIaYiIiNwLi5GHiQ7yRa9Qf5gtAvJOXhA7DhERkVthMfJAoxNb5xnxcBoREdG1YDHyQKM5z4iIiKhLWIw80IjeIZBJJThdpce5mkax4xAREbkNFiMPpFb5IDkmCAAPpxEREV0LFiMPNcp6OO0nHk4jIiLqNBYjDzXGOgF7+4kqWCyCyGmIiIjcA4uRh0qOCYKfQoaaRhOKyuvFjkNEROQWWIw8lI9MitS4YADAdq5nRERE1CluX4yWLVuGuLg4qFQqpKenY9euXZ163McffwyJRIIpU6Y4NqCIMnqFAAAXeiQiIuokty5Ga9euRXZ2NhYsWICCggIkJSVhwoQJqKio6PBxZ86cwZ///GeMGTPGSUnFMbJ3SzHaefoCzJxnREREdFVuXYxef/11zJgxA1lZWRgwYABWrFgBPz8/rF69+oqPMZvNePDBB/H888+jV69eHT6/wWCATqdrc3MnA6PUCFTKUd/UjMNl7pWdiIhIDG5bjIxGI/Lz85GZmWnbJpVKkZmZiby8vCs+buHChQgPD8cjjzxy1dfIycmBRqOx3WJiYuyS3VnkMinS4lvmGeWd4mn7REREV+O2xaiqqgpmsxkRERFttkdERECr1bb7mG3btuHf//43Vq1a1anXmDt3Lurq6my3kpKS687tbBm9Oc+IiIios+RiB3CW+vp6PPTQQ1i1ahVCQ0M79RilUgmlUungZI7VWox2na6GyWyBj8xtuzAREZHDuW0xCg0NhUwmQ3l5eZvt5eXliIyMvGz/kydP4syZM5g8ebJtm8ViAQDI5XIUFRWhd+/ejg0tgv6RagT5+aC20YQDpXUYFttN7EhEREQuy22HDxQKBVJSUpCbm2vbZrFYkJubi4yMjMv279evHw4cOIDCwkLb7Y477sC4ceNQWFjodvOHOksqlSC9dZ4RD6cRERF1yG1HjAAgOzsb06dPR2pqKtLS0rB06VLo9XpkZWUBAKZNm4bo6Gjk5ORApVJh0KBBbR4fFBQEAJdt9zQZvUKw8VA5dpy6gMfHJYgdh4iIyGW5dTGaOnUqKisrMX/+fGi1WiQnJ2PDhg22Cdlnz56FVOq2g2J2k9G7ZU7VnjM1MDZboJDz74SIiKg9EkEQuPJfJ+l0Omg0GtTV1UGtVosdp9MsFgEpL25CTaMJnz82kvOMiIjIq1zL5zeHDryAVCrBcOt103afrhY5DRERketiMfISrQs97mIxIiIiuiIWIy+RHm9dz+hMNa+bRkREdAUsRl6if/dA+CtkqG9qRpG2Xuw4RERELonFyEvIZVKkxLUeTuN6RkRERO1hMfIirQs97j5TI3ISIiIi18Ri5EVaJ2DvPF0NrtJARER0ORYjLzKkhwYKuRRVDQacrtKLHYeIiMjlsBh5EaVchuSYIAA8bZ+IiKg9LEZeJp3rGREREV0Ri5GX+eU8IyIiImqLxcjLDIvtBplUgtLaiyitvSh2HCIiIpfCYuRl/JVyDIpquYAer5tGRETUFouRF+LhNCIiovaxGHmhtNbrpnEFbCIiojZYjLzQ8LhuAICTlXpUNRhETkNEROQ6WIy8UJCfAv0iAwFwnhEREdEvsRh5qeGtF5Q9w2JERETUisXIS6VxoUciIqLLsBh5qdZidPi8Dromk8hpiIiIXAOLkZeKUKsQF+IHQQDyz9SIHYeIiMglsBh5Ma5nRERE1BaLkRdrnYC9mxOwiYiIALAYebV060KP+8/V4qLRLHIaIiIi8bEYebGYYF+EByphMgvYd65W7DhERESiYzHyYhKJBMOt84z28HAaERERi5G3G96z5fIgu3lmGhEREYuRt0u1TsAuKK6B2SKInIaIiEhcLEZerl9kIAKUctQbmlGkrRc7DhERkahYjLycXCbF0NggAMCeYs4zIiIi78ZiRL9Yz4jzjIiIyLuxGBFS46wTsE9XQxA4z4iIiLwXixEhOSYIcqkEWl0TSmsvih2HiIhINCxGBD+FHAOjNQCAPTycRkREXozFiAD8cj0jTsAmIiLvxWJEAC6tZ8QRIyIi8mYsRgTg0gTsovJ61DYaRU5DREQkDhYjAgCEBijRK8wfAJBfzFEjIiLyTixGZDO8J9czIiIi78ZiRDath9P2cAI2ERF5KRYjsmldAXv/uTo0mcwipyEiInI+FiOy6Rnih9AAJYxmCw6U1okdh4iIyOlYjMhGIpFgeBzXMyIiIu/FYkRtcD0jIiLyZm5fjJYtW4a4uDioVCqkp6dj165dV9z3888/R2pqKoKCguDv74/k5GR88MEHTkzr+ob/YgK2xcILyhIRkXdx62K0du1aZGdnY8GCBSgoKEBSUhImTJiAioqKdvcPDg7GX//6V+Tl5WH//v3IyspCVlYWNm7c6OTkrmtAdzX8FDLomppxvKJB7DhEREROJREEwW2HBdLT0zF8+HC8/fbbAACLxYKYmBg88cQTePbZZzv1HMOGDcOkSZPwwgsvXHafwWCAwWCwfa/T6RATE4O6ujqo1Wr7vAkX9OC/duDnExfw4pRB+N2InmLHISIiui46nQ4ajaZTn99uO2JkNBqRn5+PzMxM2zapVIrMzEzk5eVd9fGCICA3NxdFRUW44YYb2t0nJycHGo3GdouJibFbfleW2rN1nhEnYBMRkXdx22JUVVUFs9mMiIiINtsjIiKg1Wqv+Li6ujoEBARAoVBg0qRJeOutt3DzzTe3u+/cuXNRV1dnu5WUlNj1PbiqtHiugE1ERN5JLnYAZwsMDERhYSEaGhqQm5uL7Oxs9OrVCzfeeONl+yqVSiiVSueHFFlyTBBkUglKay+irPYiooJ8xY5ERETkFG5bjEJDQyGTyVBeXt5me3l5OSIjI6/4OKlUioSEBABAcnIyjhw5gpycnHaLkbfyV8oxMEqN/efqsPtMNe5MjhY7EhERkVO47aE0hUKBlJQU5Obm2rZZLBbk5uYiIyOj089jsVjaTLCmFpfmGfFwGhEReQ+3HTECgOzsbEyfPh2pqalIS0vD0qVLodfrkZWVBQCYNm0aoqOjkZOTA6BlMnVqaip69+4Ng8GAb775Bh988AGWL18u5ttwScPjumH1z6e5AjYREXkVty5GU6dORWVlJebPnw+tVovk5GRs2LDBNiH77NmzkEovDYrp9Xo89thjOHfuHHx9fdGvXz98+OGHmDp1qlhvwWWlWBd6LCqvR91FEzS+PiInIiIicjy3XsfI2a5lHQRPcOOSzThzoRHvZg3HuL7hYschIiLqEq9Yx4gc79J103g4jYiIvAOLEV1R63XTuJ4RERF5CxYjuqLWEaN9JbUwNJtFTkNEROR4LEZ0Rb1C/RHsr4Ch2YKDpTqx4xARETkcixFdkUQiQWrPlsNpnGdERETegMWIOsTrphERkTdhMaIOtc4zyi+uhsXClR2IiMizsRhRhwZGqaHykaKm0YSTlQ1ixyEiInIoFiPqkI9MiqExPG2fiIi8A4sRXVXrekacgE1ERJ6OxYiuqnWe0e5iFiMiIvJsLEZ0VUNjgyCVACXVF6GtaxI7DhERkcOwGNFVBap80L97y0X39nDUiIiIPBiLEXXKcNsFZTkBm4iIPBeLEXVKqu2CshwxIiIiz8ViRJ3SOmJ05LwO9U0mkdMQERE5BosRdUqEWoXYYD9YBGDv2Vqx4xARETkEixF1WirXMyIiIg/HYkSd1no4jStgExGRp2Ixok5rXQF7b0kNTGaLyGmIiIjsj8WIOq13WAC6+fmgyWTBwdI6seMQERHZHYsRdZpEIkFKT65nREREnovFiK7JcK5nREREHozFiK5J6wVl9xTXQBAEkdMQERHZF4sRXZNB0Woo5VJU6404VaUXOw4REZFdsRjRNVHKZUiKCQLA9YyIiMjzsBjRNbs0z4gTsImIyLOwGNE1s80z4ogRERF5GBYjumYpPbtBIgHOXGhERX2T2HGIiIjshsWIrpla5YN+kWoAQD4PpxERkQdhMaIu4TwjIiLyRCxG1CWX1jPiPCMiIvIcLEbUJa0jRofKdNAbmkVOQ0REZB8sRtQl3TW+iA7yhdkiYO/ZWrHjEBER2QWLEXUZr5tGRESehsWIuozzjIiIyNOwGFGXDbcWo71na2EyW0ROQ0REdP1YjKjLEsMDoFbJ0Wg048h5ndhxiIiIrhuLEXWZVCqxHU7jekZEROQJWIzougznddOIiMiDXHcxqq6uhsXC+SXe6pcrYAuCIHIaIiKi69OlYnT48GEsWrQII0eORFhYGMLDwzFt2jR89tln0Ov19s5ILmxwDw0UcimqGgwovtAodhwiIqLr0uliVFRUhKeffhqJiYkYMWIEdu/ejUcffRTl5eX45ptv0LNnTyxcuBChoaGYOHEili9f7sjc5CKUchmSemgAcD0jIiJyf50uRtu3b4der8ebb76JqqoqfPbZZ5g2bRpCQ0ORlpaGF154Afv27cORI0dw66234rPPPnNkbnIhtvWMOAGbiIjcXKeLUVZWFlasWIGJEydCoVC0ua+goMD257i4OMyaNQvff/+9/VJ2YNmyZYiLi4NKpUJ6ejp27dp1xX1XrVqFMWPGoFu3bujWrRsyMzM73J86hytgExGRp7DLWWlpaWnIzs5us+2bb76xx1N3aO3atcjOzsaCBQtQUFCApKQkTJgwARUVFe3uv2XLFjzwwAPYvHkz8vLyEBMTg1tuuQWlpaUOz+rJUmKDIZEAp6r0qKw3iB2HiIioy+xSjAYPHgy1Wo2srCzbtr/97W/2eOoOvf7665gxYwaysrIwYMAArFixAn5+fli9enW7+3/00Ud47LHHkJycjH79+uFf//oXLBYLcnNzHZ7Vk2n8fNA/Ug0A2Hn6gshpiIiIus4uxUgikeDvf/87kpKScM8998BkMjn81G2j0Yj8/HxkZmbatkmlUmRmZiIvL69Tz9HY2AiTyYTg4OB27zcYDNDpdG1u1L70Xi1/hztP8XAaERG5L7sUI7W6ZbRg9uzZmDx5Mu644w5cvHjRHk99RVVVVTCbzYiIiGizPSIiAlqttlPPMWfOHERFRbUpV7+Uk5MDjUZju8XExFx3bk81olcIAGDHKY4YERGR++pSMWpoaGjz/ZYtW2x/nj59Ov74xz9ecZ6Pq1i0aBE+/vhjrFu3DiqVqt195s6di7q6OtutpKTEySndR5r1zLTjFQ2oauA8IyIick9dKkYajabD0/HvuusuVFc79pBKaGgoZDIZysvL22wvLy9HZGRkh4999dVXsWjRInz33XcYMmTIFfdTKpVQq9VtbtS+bv4K9IsMBADsOs3DaURE5J66VIwEQcDKlSsxatQojB49GrNnz8bu3bvtna1DCoUCKSkpbSZOt06kzsjIuOLjXnnlFbzwwgvYsGEDUlNTnRHVa/BwGhERubsuzzHau3cvhg0bhtGjR+PQoUMYM2YM/vznP9sz21VlZ2dj1apVeP/993HkyBH86U9/gl6vt50dN23aNMydO9e2/+LFizFv3jysXr0acXFx0Gq10Gq1lx0apK4ZwQnYRETk5uRdfeCaNWtw8803277fv38/7rzzTkRHR+Opp56yS7irmTp1KiorKzF//nxotVokJydjw4YNtgnZZ8+ehVR6qfstX74cRqMR99xzT5vnWbBgAf7+9787JbMnS4tvGTEqKq9Htd6IYH/FVR5BRETkWiRCF86rDw0NxbZt29CvX7822//3v//hqaeewrFjx+wW0JXodDpoNBrU1dVxvtEVTPjHjygqr8eK3w3DrYO6ix2HiIjomj6/u3QoLTk5Ge++++5l2xMSEnD27NmuPCV5iNb1jHbwcBoREbmhLhWjF198EW+++SYeeugh5OXlQa/Xo6KiAi+//DLi4+PtnZHcCCdgExGRO+vSHKMRI0Zgx44dmDVrFsaMGWNb5VqlUuGTTz6xa0ByL2nxLSNGR7X1qNEb0Y3zjIiIyI10efJ1UlIStmzZgoqKCuTn58NisSA9PR2hoaH2zEduJjRAicTwAByvaMCuM9WYMLDjNaWIiIhcSacPpWm1WhgMl69oHB4ejokTJ2LSpEltStGpU6fsk5DczqV5RjycRkRE7qXTxejTTz9FcHAw7rrrLrz77ruorKy8bJ+dO3fiueeew8CBA5GUlGTXoOQ+WucZcT0jIiJyN50uRjNnzsS+ffswZswYvPfee+jRowdGjx6Nl19+GTNmzED37t0xZcoUVFRUYNGiRe0WJ/IOrfOMjmh1qGs0iZyGiIio87q0jhEAXLhwAevXr8c333yDuLg43HnnncjIyIBEIrF3RpfBdYw6b/xrW3CyUo9V01Jx84AIseMQEZEXu5bP7y5Pvg4JCcH06dMxffr0rj4FebD0XiE4WanHjlMXWIyIiMhtdPlaaUQdsc0zOs0J2ERE5D5YjMghRljnGR0q06HuIucZERGRe2AxIocIV6vQK9QfggDsOcOz04iIyD2wGJHDcD0jIiJyNyxG5DCX5hlxxIiIiNwDixE5THp8SzE6WFrHeUZEROQWWIzIYSI1KvQK84dFAHbycBoREbkBFiNyqFG9W66f9/OJKpGTEBERXR2LETnUqARrMTrJESMiInJ9LEbkUBm9QiCVACcqGqCtaxI7DhERUYdYjMihNH4+GBytAQBsP8nDaURE5NpYjMjhRloPp23jPCMiInJxLEbkcK0TsLefuABBEEROQ0REdGUsRuRwqXHdoJBLodU14WSlXuw4REREV8RiRA6n8pEhtWc3AJxnREREro3FiJzCdto+5xkREZELYzEip2gtRnknL8Bs4TwjIiJyTSxG5BSDozUIVMmha2rGwdI6seMQERG1i8WInEImlSCjV8tFZXnaPhERuSoWI3Ka1sNpnIBNRESuisWInKa1GO0+U4Mmk1nkNERERJdjMSKn6R3mjwi1EsZmC/KLa8SOQ0REdBkWI3IaiURiGzXiPCMiInJFLEbkVJcuD8JiRERErofFiJyqdcRof2kd6hpNIqchIiJqi8WInCpSo0LvMH8IApB36oLYcYiIiNpgMSKn4+VBiIjIVbEYkdPZihHXMyIiIhfDYkRON6JXCKQS4FSlHqW1F8WOQ0REZMNiRE6n8fVBckwQAODHY5XihiEiIvoFFiMSxdg+4QBYjIiIyLWwGJEobuhzaaHHZrNF5DREREQtWIxIFEN6BCHIzwf1Tc3Yd65W7DhEREQAWIxIJDLppcuDbC3i4TQiInINLEYkmrGJYQCArcd52j4REbkGty5Gy5YtQ1xcHFQqFdLT07Fr164r7nvo0CH85je/QVxcHCQSCZYuXeq8oNSuG/q0FKP952pRozeKnIaIiMiNi9HatWuRnZ2NBQsWoKCgAElJSZgwYQIqKira3b+xsRG9evXCokWLEBkZ6eS01J5IjQp9IwIhCC2TsImIiMTmtsXo9ddfx4wZM5CVlYUBAwZgxYoV8PPzw+rVq9vdf/jw4ViyZAnuv/9+KJVKJ6elK2k9O20rT9snIiIX4JbFyGg0Ij8/H5mZmbZtUqkUmZmZyMvLs9vrGAwG6HS6Njeyr9bDaT8eq4QgCCKnISIib+eWxaiqqgpmsxkRERFttkdERECr1drtdXJycqDRaGy3mJgYuz03tRgeFwxfHxkq6g04cr5e7DhEROTl3LIYOcvcuXNRV1dnu5WUlIgdyeOofGQYlRACANhc1P78MCIiImdxy2IUGhoKmUyG8vLyNtvLy8vtOrFaqVRCrVa3uZH93di35fIgm4+yGBERkbjcshgpFAqkpKQgNzfXts1isSA3NxcZGRkiJqOuGNevpRgVnK3haftERCQqtyxGAJCdnY1Vq1bh/fffx5EjR/CnP/0Jer0eWVlZAIBp06Zh7ty5tv2NRiMKCwtRWFgIo9GI0tJSFBYW4sSJE2K9BbKKDvJF34hAWATgx+M8O42IiMQjFztAV02dOhWVlZWYP38+tFotkpOTsWHDBtuE7LNnz0IqvdT7ysrKMHToUNv3r776Kl599VWMHTsWW7ZscXZ8+pVx/cJRVF6PzUcrcGdytNhxiIjIS0kEniPdaTqdDhqNBnV1dZxvZGc7T13A1Hd2oJufD/b87WbIpBKxIxERkYe4ls9vtz2URp4lpWc3BKrkqGk0obCkVuw4RETkpViMyCXIZVLbYo9beNo+ERGJhMWIXMZN1tP2f+Bp+0REJBIWI3IZY/uGQSIBDpXpoK1rEjsOERF5IRYjchmhAUoMjQkCAOQeLe94ZyIiIgdgMSKXkjmgZbmF7w+zGBERkfOxGJFLubl/SzH6+eQF6A3NIqchIiJvw2JELiUhPAA9Q/xgbLbgp+NVYschIiIvw2JELkUikSDTOmr0/REeTiMiIudiMSKX01qMfjhaAbOFC7MTEZHzsBiRy0mN6waNrw+q9UbsPVsjdhwiIvIiLEbkcnxkUozr27IK9iYeTiMiIidiMSKX1Hra/iaetk9ERE7EYkQu6YY+YfCRSXCqUo8TFQ1ixyEiIi/BYkQuSa3ywaiEUADAxkNakdMQEZG3YDEil3XrwEgAwIaDLEZEROQcLEbksjIHREAqAQ6U1qGkulHsOERE5AVYjMhlhQYoMTwuGAAPpxERkXOwGJFLu3VQy+E0FiMiInIGFiNyaROs84z2FNegor5J5DREROTpWIzIpUUF+SIpJgiCAHx3iGsaERGRY7EYkctrPTuNh9OIiMjRWIzI5bXOM9p+8gJq9EaR0xARkSdjMSKXFx/qj/7d1TBbBGzgqBERETkQixG5hduHdAcArN9fJnISIiLyZCxG5BYmD4kCAOSdvIDKeoPIaYiIyFOxGJFbiA3xw5AeGlgEYMPB82LHISIiD8ViRG6j9XDa1/tZjIiIyDFYjMhtTLIeTtt9phrlOi72SERE9sdiRG4jOsgXw2JbFnv8H0eNiIjIAViMyK3cbh014tlpRETkCCxG5FYmDekOiQQoOFuLkupGseMQEZGHYTEitxKhViGjVwgA4MvCUpHTEBGRp2ExIrdz19BoAMDne0shCILIaYiIyJOwGJHbuXVQJJRyKU5V6nGgtE7sOERE5EFYjMjtBKp8cMvAlgvLfl7Aw2lERGQ/LEbklu4a2nJ22tf7ymAyW0ROQ0REnoLFiNzSmMQwhPgrcEFvxLYTVWLHISIiDyEXOwBRV/jIpJicFIX3tp/BuoJSjOsbLnYkIo8kCALqLppQWW/ABb0Ruosm6JqaUd9kgu5iM3RNJjQ0NaOp2QyDyQJDsxmGZguaTGYYfzGaK4Gk5asEv9jWskEpl0LlI4NKLoXS+lXlI4PKp/WrDAFKectNJUegSo5ApQ8CVC3bAlVyKOVSSH755ERdxGJEbuuuodF4b/sZfHdYi/omEwJVPmJHInIrgiDggt6IczUXca6m0fa1rLYJlfUGaxkywGR2/bM/FTIpgvx8EOyvQDc/BYL9W27d/BUI9vNBN38FQvyV6OZ/aR+Vj0zs2OSCWIzIbQ3poUFCeABOVDTg633n8dv0WLEjEbkkY7MFxRf0OFHR0HKrbMDx8gacrtLjosncqedQq+QIDVBC7evTclPJEajygdpXjkClHCofGZQ+MijlUuut5c9oHcT5RbcSfvGNxQLbCFPrqFNTsxlNJgsMppbRp4tGMxoMzag3NKOhyYT6pmY0GJrR0NSMBmMzBAEwmi2oqDegot7Q6b+XQJUc4YFKhAUqERaoQliAEuFqJcICWrcpER6oRDc/BaRSjkZ5CxYjclsSiQRTU2Pw0jdHsHZPCYsREYAmkxlHzutwsLQOB0rrcKBUh+Pl9Wi2tD/qI5EAEYEq9Ojma735ISrIFxFqJUIDlAgNVCI0QAGl3DVHVywWAXpjM3RNzajRG1HTaES1vuVWozfiwq+2VetNqGk0wmwRUN/UjPqmZpys1Hf4GjKpBKEBCmtRUtmKU2uJCldbtwcqOQrlAViMyK3dNSwaizccxb6SWhRp69E3MlDsSEROVVZ7EbvPVGPX6WrkF9fgeEUDzO2UIH+FDAnhAegdHoDE8MCWP4f5I7qbr8uWns6QSiUIVPkgUOWD6CDfTj2mdd5UVUPLCFPrYcPKBgMqddavtkOJLSWqXGdAuc4AQNfhc2t8fRAeeKks2b4GKhGhVtnu81Pw49dV8SdDbi00QInx/cOx8VA51u4uwfzJA8SORORQJdWN+PlEFXZZy9C5mouX7RMaoMCgaA0GRWlavkarER3ky8nJVhKJBEF+CgT5KZAQ3vH/TJnMFlxoMFqL06W5VxW/+lqua4Kh2YK6iybUXTTheEVDh88boJRbS1NLcYr4dZGy3heglPPn5mRuX4yWLVuGJUuWQKvVIikpCW+99RbS0tKuuP8nn3yCefPm4cyZM0hMTMTixYtx2223OTEx2dvU4THYeKgc6/aew5yJfd36/36Jfq3JZMau09XYUlSJrccqLjvsI5UAg6I1GB4XjOFx3ZAUE4RItYofpnbiI5MiUqNCpEYFQHPF/QRBgK6pGRW6JutcpyaU6wyo0LX8ufVruc6Ai6aWOVMNlc04dZXDeL4+MltpClMrEeqvQLC/EsEBCgRbJ5mHBLR8DfL1gVzGVXiul1sXo7Vr1yI7OxsrVqxAeno6li5digkTJqCoqAjh4Zefvr19+3Y88MADyMnJwe233441a9ZgypQpKCgowKBBg0R4B2QPNySGIUKtRLnOgO8PV2DSkO5iRyK6LpX1Bnx3WItNh8ux49QFNJkunfYuk0owNCYII3qFIC0+GMN6dkOA0q1/lXsEiUQCja8PNL4+SIy48iiUIAhoMDS3lKd2StOl7w1oMDTjosmMMxcaceZCYycytBzKC/ZXIMR65l1raQr2VyLY3wfB/kqE+CugsU6iD1TKObH8VySCG1+FMz09HcOHD8fbb78NALBYLIiJicETTzyBZ5999rL9p06dCr1ej/Xr19u2jRgxAsnJyVixYsVVX0+n00Gj0aCurg5qtdp+b4Su25KNR7Fs80nc0CcM/3n4yiOGRK5KW9eEDQfP45uDWuw+U41f/maOVKswtk8YbuwbhpEJodD4cmkKb6C3FaiWUahyXVPLpPJGIy40WCeTWyeW1zaauvQaEknLYT2Nrw/U1rMML/3Zx/pn+aU/+/rAXyGHv1IGP+tXXx+Zy49QXsvnt9v+b4bRaER+fj7mzp1r2yaVSpGZmYm8vLx2H5OXl4fs7Ow22yZMmIAvvvii3f0NBgMMhkunfup0HU+6I/HcmxKDZZtP4qfjlSipbkRMsJ/YkYiuqlpvxFeFpfhqXxkKzta2uS+phwYTBkVifL8I9IkIcPkPHrI/f6Uc8Uo54kP9r7pvs9mCmkbT5aWpwYhqvcF2dl7rfbomE5pMFggCbGfnAZfPV+sMiQTw85HBTymHv+IXhUnR9ns/hfzSop3yS4t3qnyktqUeVD4yqFU+SAgP6FIWe3DbYlRVVQWz2YyIiIg22yMiInD06NF2H6PVatvdX6vVtrt/Tk4Onn/+efsEJoeKC/XHqIQQ/HziAtbsOos5t/YTOxJRu4zNFmwuqsBn+eewuaiizeKJKT27YeKgSNw6KBI9urHcU+fJZVLb2kuIuPr+AGBoNttWL9dZJ43rmpp/8WdTm/tbt+uNZjQamtFoMkMQAEEA9EYz9EYzKu3wXvp3V+PbWWPs8Exd47bFyBnmzp3bZoRJp9MhJiZGxETUkYdGxOHnExewdncJZmcmchI2uZQj53VYu7sEX+0rQ7XeaNs+KFqNu4f2wG2Du1sn+BI5h1IuQ1igrKVMdYHFIqCp2Qy9wYxGYzMajS1fW7+3fbUWqQaD2bp4p3UhT5O5zeKeTdZt3UX+d+C2xSg0NBQymQzl5eVttpeXlyMyMrLdx0RGRl7T/kqlEkpl1/6DIefL7B+OSLUKWl0Tvj2gxZSh0WJHIi9naDZjw0EtPsgrxp7iGtv28EAl7hoajbuH9eDaW+S2pFIJ/BRy65pMnvNZ6bbn9SkUCqSkpCA3N9e2zWKxIDc3FxkZGe0+JiMjo83+ALBp06Yr7k/uRS6T2la//mBHschpyJudq2nEKxuOYmTOD5j1cSH2FNdALpVg0uDueC9rOLY/exPm3tafpYjIBbntiBEAZGdnY/r06UhNTUVaWhqWLl0KvV6PrKwsAMC0adMQHR2NnJwcAMCsWbMwduxYvPbaa5g0aRI+/vhj7NmzB++8846Yb4Ps6P60GLyZexz5xTU4VFaHgVFXXneEyN7yi2vwzo8nselwOVoXn45Uq/Db9FjcPzwG4WoeKiNydW5djKZOnYrKykrMnz8fWq0WycnJ2LBhg22C9dmzZyGVXhoUGzlyJNasWYO//e1veO6555CYmIgvvviCaxh5kPBAFW4dFIn1+8/jwx1nkXP3YLEjkYezWATkHq3Ayq0n2xwuG50Qit+N6InM/uFcdI/Ijbj1OkbOxnWM3MPOUxcw9Z0d8PWRYcfc8dD4cc0Xsj9DsxnrCkqx6qdTttWofWQSTEmOxowbeqFPB4v8EZFzecU6RkRXkhYfjP7d1ThyXoePdhXjsRsTxI5EHqTJZMZ/95Tgn5tPQqtrAgAEquR4ML0nskbFIYKHy4jcGosReRyJRII/jI7H05/sw/vbz+APo3tBIeehDLo+TSYz1u4uwfItlwpRpFqFP4yJx/1psbwsB5GH4L9k8kiTk6KweMNRlOsM+HpfGX6T0kPsSOSmWgvRP7ecQLmuZSX87hoVHhuXgPtSe3C9LCIPw2JEHkkhl2L6yDgs2ViEf207jbuHRfOSCnRNzBYBnxecwz82HUNZXcsIEQsRkedjMSKP9WB6LN7+4QSOnNdh+8kLGJUQKnYkcgOCIOCHoxVYvOEojpU3AGgpRI+PS8C9LEREHo/FiDxWkJ8C96X2wPt5xVj10ykWI7qq/OJqLPr2KHafaTntXuPrg8fH9ca0jDiofFiIiLwBixF5tIdHx+ODHcXYUlSJw2U6DIjiMgt0ubMXGpHz7RF8e7DlgtJKuRRZo+Lxp7G9udwDkZdhMSKP1jPEH7cN7o71+8/j7c3H8c8HU8SORC6kwdCMZZtP4N8/nYbRbIFUAtyXGoNZmYnorvEVOx4RiYDFiDzezJsSsH7/eXx7UIvj5fVI5MJ7Xs9iEfBpwTks2ViEyvqWM81GJ4Ri3u0DeP0yIi/HYkQer1+kGhMGRmDjoXK8vfkE3rh/qNiRSET5xdX4+1eHcaC0DgAQF+KHv04agMz+4TxzkYhYjMg7PHFTIjYeKsfX+8owO7MP4kP9xY5ETlajN2LRt0exdk8JACBQKccT4xMwfWQczzQjIhsuB0xeYVC0Bjf1C4dFAJZtPiF2HHIii0XA2t1ncdNrW2yl6N6UHtj8zI344w29WYqIqA2OGJHXeOKmBPxwtALr9pbi8XEJHDXyAkfO6/C3Lw4i33rV+36RgXhxyiCkxgWLnIyIXBVHjMhrDI3thnF9w2C2CHjtuyKx45ADNRia8eL6w7j9rW3IL66Bn0KGv03qj6+fGM1SREQd4ogReZVnJvTDlmOVWL//PP7fDXUY3EMjdiSyI0EQ8O1BLRZ+fdh2odfbBkdi3u0DePo9EXUKR4zIqwyIUuPOpCgAwCsbj4qchuzpfN1FPPL+Hjz2UQG0uibEBvvh3azh+OeDKSxFRNRpHDEir/P0LX3xvwPn8dPxKvx8ooqXCnFzgiDg/3aVIOebI6g3NEMhk+LRsb3w2LgEXsaDiK4ZR4zI68QE++HB9J4AgMUbjsJiEURORF119kIjHvzXTjy37gDqDc0YGhuEb2aNRvYtfVmKiKhLWIzIK828KQH+Chn2n6vD53tLxY5D18hsEbB622lMWPojtp+8AJWPFH+b1B+fPjoSCeFcuZqIuo7FiLxSaIAST4xPBAAs+vYo6ptMIieizjpR0YD7VuZh4frDuGgyY0SvYGyYdQP+MKYXZFKuXE1E14fFiLxW1qg4xIf6o6rBgDdzj4sdh66i2WzBP7ecwG1v/oT84hoEKOV4ccogrPnDCMRxTSoishMWI/JaSrkM8ycPAAC8+/MZnKioFzkRXcmR8zrc9c/teGVDEYzNFoztE4aNT92A343oCSlHiYjIjliMyKuN6xuOzP7haLYIeP7rwxAETsR2JcZmC17fdAyT39qGA6V1UKvkePXeJLyXNRzRQTwFn4jsj8WIvN682wdAIZPip+NV+GpfmdhxyGpfSS0mv7UNb+YeR7NFwC0DIvB99ljck9IDEglHiYjIMViMyOv1DPHHEzclAAD+/tUhXGgwiJzIuzWZzMj55gju+ufPKCqvR4i/Am//dihWPpSCcLVK7HhE5OFYjIgA/L+xvdEvMhA1jSY8//VhseN4rd1nqjHxjZ+w8sdTsAjAnclR2JQ9FrcPieIoERE5BYsREQCFXIpX7hkCqQT4al8Zco+Uix3Jq+gNzVjw5UHctzIPp6v0CA9UYtW0VLxx/1AE+yvEjkdEXoTFiMhqSI8gzBjTCwDw13UHUdfItY2cYdvxKkxY+iPezyuGIABTU2OwKXssbh4QIXY0IvJCLEZEvzA7sw/iQ/2h1TXh2c/38yw1B6q7aMKcT/fjd//eiXM1FxEd5IsPHknD4nuGQOPrI3Y8IvJSLEZEv+CrkOGN+5PhI5Pg24NafLy7ROxIHmnT4XLc8o+tWLun5e93WkZPbHzqBoxJDBM5GRF5OxYjol8Z0iMIf76lLwDg+a8PceFHO6pqMGDmmgLM+M8elOsMiAvxw9o/jsDCOwchQCkXOx4REYsRUXtmjOmFMYmhaDJZMHPNXlw0msWO5NYEQcAXe0tx8+tbsX7/ecikEjw6tjc2zL4B6b1CxI5HRGTDYkTUDqlUgtfuTUKIvwJHtfWY8xnnG3VVWe1FPPzebsxeW4iaRhP6d1fji8dG4dmJ/aDykYkdj4ioDRYjoisIV6vw9m+HQSaV4Kt9ZXjnx1NiR3IrFouAD3YU45Z//IjNRZVQyKT48y198NXMURjcQyN2PCKidrEYEXUgo3cIFlgvNLt4w1FsPVYpciL3cKqyAfev2oF5XxxEg6EZw2KD8M2s0Zh5UyJ8ZPy1Q0Sui7+hiK7ioRE9MTU1BhYBeGJNAY6VczL2lTSZzFj6/THc+sZP2HW6Gr4+MiyYPACfPDoSCeGBYscjIroqngZCdBUSiQQLpwzEycoG7CmuwfTVu/DZn0Yiild3b2Pb8SrM+/IgTlfpAQBjEkPx8l2DERPsJ3IyIqLO44gRUSco5TL8a3oqEsIDcL6uCdNX70Jto1HsWC6hor4Jsz7ei9/9e6ftch5v/3Yo/vNwGksREbkdFiOiTgryU+D9h9MQoVbieEUDHnl/DxoMzWLHEk2z2YL3t5/B+Ne24svCMkgkwO9HxuH7p3nRVyJyXxKB5yB3mk6ng0ajQV1dHdRqtdhxSCRHtTrcuyIP9U0tk4rfzUrzuktY/HisEi+sP4zjFQ0AgMHRGrx01yAM6REkbjAionZcy+c3i9E1YDGiVvtKajFt9S7UXTRhULQaHzycjm5ecBX4U5UNeOl/R5B7tAIA0M3PB9k398Fv03tCJuUIERG5JhYjB2Exol86XKbDQ//eiQt6I/pGBOLdrOEeOyFb12TCW7nH8d72MzCZBcilEjyU0ROzx/eBxs+7RsuIyP2wGDkIixH92omKevx21U5U1BsQFqjEqmmpSI4JEjuW3TQam/He9jNYufUU6i6aAAA39g3D3yYNQEJ4gMjpiIg6h8XIQViMqD2ltRfxyHu7cVRbD6VcitfuS8LtQ6LEjnVdmkxmfLTzLJZvOYGqhpaz7xLCA/DX2/pjXL9wkdMREV0bFiMHYTGiK2kwNGPW/+21zb2ZltETz93W3+2uBWZstmDtnhIs++EEtLomAEBssB9mZybizuRoziMiIrd0LZ/fbnm6fnV1NR588EGo1WoEBQXhkUceQUNDQ4ePeeedd3DjjTdCrVZDIpGgtrbWOWHJKwQo5XhnWir+39heAID/5BXjzrd/dptVsmv0RizbfAJjXvkB8744CK2uCVEaFXLuHozcp8fi7mE9WIqIyCu45YjRxIkTcf78eaxcuRImkwlZWVkYPnw41qxZc8XHLF26FE1NLf8HPHfuXNTU1CAoKOiaXpcjRtQZW49V4un/FqKqwQiFXIo/je2NP93Y2yVHj05U1GP1z2fwecE5NJksAICwQCVmjkvA/WkxUMpdLzMR0bXy6ENpR44cwYABA7B7926kpqYCADZs2IDbbrsN586dQ1RUx3M7tmzZgnHjxnWqGBkMBhgMBtv3Op0OMTExLEZ0VZX1Bvzl033YXNRy0dmeIX5YMHkAxvUNF33hw2azBT8er8T724vbXBR3QHc1Hhkdj9uTurMQEZFHuZZi5HbXSsvLy0NQUJCtFAFAZmYmpFIpdu7cibvuustur5WTk4Pnn3/ebs9H3iMsUInVvx+Obw5o8cL6wyi+0IiH39uDYbFBeHJ8Isb2CXNqQRIEAUe19fhibyk+31uKyvqWwi+RADf3j8Ajo+ORFh8semkjIhKb2xUjrVaL8PC2Z8XI5XIEBwdDq9Xa9bXmzp2L7Oxs2/etI0ZEnSGRSDBpSHeM7RtmWwOo4Gwtfv/ubgyO1uCBtFhMTuqOQJVj1gEyWwQcLK3Dd4e1+PaAFqesF3cFgGB/Be4aGo3pGXGIDeH1zIiIWrlMMXr22WexePHiDvc5cuSIk9K0UCqVUCqVTn1N8jwBSjnm3tYfj4yJx8qtp/DRzmIcKK3DgXUHsHD9IdwyIBI39g3D6MRQhAequvw6JrMFRdp67C2pxY5TF/DziSrUNpps9yvkUtzYJwy/SemBcX3DoZC75bkXREQO5TLF6Omnn8bvf//7Dvfp1asXIiMjUVFR0WZ7c3MzqqurERkZ6cCERNcnPFCFebcPwJ9u7I3PC85h7e4SnKzU46t9ZfhqXxkAoHeYP/p3V6NfZCDiQwMQ7K9AN38f+CvkEATAIghoNJpxQW/AhQYjSmsv4mRFA05WNqCovN42gbpVoFKOUQmhmDg4EuP7RyBA6TL/5ImIXJLL/JYMCwtDWFjYVffLyMhAbW0t8vPzkZKSAgD44YcfYLFYkJ6e7uiYRNctNECJP97QGzPG9MLeklp8f7gcPx2vwsGyOpys1ONkpR7r95/v0nMHquRIjgnCsNhuuKFPKJJ6BEEu48gQEVFnuUwx6qz+/fvj1ltvxYwZM7BixQqYTCbMnDkT999/v+2MtNLSUowfPx7/+c9/kJaWBqBlbpJWq8WJEycAAAcOHEBgYCBiY2MRHBws2vsh7yWRSDAsthuGxXbDX24FqvVG7DtXiyJtPY6e1+FczUXUNBpR02hCo7EZMokEUokESh8pQgOUCAlQICJQhd7hAegV6o/EiED0CvWHlOsNERF1mdsVIwD46KOPMHPmTIwfPx5SqRS/+c1v8Oabb9ruN5lMKCoqQmNjo23bihUr2pxhdsMNNwAA3n333asewiNyhmB/Bcb1Dce4vrzkBhGRWNxuHSMxcYFHIiIi9+PxlwQhIiIicgQWIyIiIiIrFiMiIiIiKxYjIiIiIisWIyIiIiIrFiMiIiIiKxYjIiIiIisWIyIiIiIrFiMiIiIiKxYjIiIiIisWIyIiIiIrFiMiIiIiKxYjIiIiIiu52AHciSAIAFqu0ktERETuofVzu/VzvCMsRtegvr4eABATEyNyEiIiIrpW9fX10Gg0He4jETpTnwgAYLFYUFZWhsDAQEgkErs+t06nQ0xMDEpKSqBWq+363K6I79ez8f16Nm97v4D3vWdPe7+CIKC+vh5RUVGQSjueRcQRo2sglUrRo0cPh76GWq32iP8IO4vv17Px/Xo2b3u/gPe9Z096v1cbKWrFyddEREREVixGRERERFYsRi5CqVRiwYIFUCqVYkdxCr5fz8b369m87f0C3veeve39/hInXxMRERFZccSIiIiIyIrFiIiIiMiKxYiIiIjIisWIiIiIyIrFyAUsW7YMcXFxUKlUSE9Px65du8SO5DA5OTkYPnw4AgMDER4ejilTpqCoqEjsWE6xaNEiSCQSzJ49W+woDlVaWorf/e53CAkJga+vLwYPHow9e/aIHcshzGYz5s2bh/j4ePj6+qJ379544YUXOnU9Jnfw448/YvLkyYiKioJEIsEXX3zR5n5BEDB//nx0794dvr6+yMzMxPHjx8UJawcdvV+TyYQ5c+Zg8ODB8Pf3R1RUFKZNm4aysjLxAl+nq/18f+nRRx+FRCLB0qVLnZZPLCxGIlu7di2ys7OxYMECFBQUICkpCRMmTEBFRYXY0Rxi69atePzxx7Fjxw5s2rQJJpMJt9xyC/R6vdjRHGr37t1YuXIlhgwZInYUh6qpqcGoUaPg4+ODb7/9FocPH8Zrr72Gbt26iR3NIRYvXozly5fj7bffxpEjR7B48WK88soreOutt8SOZhd6vR5JSUlYtmxZu/e/8sorePPNN7FixQrs3LkT/v7+mDBhApqampyc1D46er+NjY0oKCjAvHnzUFBQgM8//xxFRUW44447REhqH1f7+bZat24dduzYgaioKCclE5lAokpLSxMef/xx2/dms1mIiooScnJyREzlPBUVFQIAYevWrWJHcZj6+nohMTFR2LRpkzB27Fhh1qxZYkdymDlz5gijR48WO4bTTJo0SXj44YfbbLv77ruFBx98UKREjgNAWLdune17i8UiREZGCkuWLLFtq62tFZRKpfB///d/IiS0r1+/3/bs2rVLACAUFxc7J5QDXen9njt3ToiOjhYOHjwo9OzZU/jHP/7h9GzOxhEjERmNRuTn5yMzM9O2TSqVIjMzE3l5eSImc566ujoAQHBwsMhJHOfxxx/HpEmT2vycPdVXX32F1NRU3HvvvQgPD8fQoUOxatUqsWM5zMiRI5Gbm4tjx44BAPbt24dt27Zh4sSJIidzvNOnT0Or1bb571qj0SA9Pd2rfn9JJBIEBQWJHcUhLBYLHnroITzzzDMYOHCg2HGchheRFVFVVRXMZjMiIiLabI+IiMDRo0dFSuU8FosFs2fPxqhRozBo0CCx4zjExx9/jIKCAuzevVvsKE5x6tQpLF++HNnZ2Xjuueewe/duPPnkk1AoFJg+fbrY8ezu2WefhU6nQ79+/SCTyWA2m/HSSy/hwQcfFDuaw2m1WgBo9/dX632erKmpCXPmzMEDDzzgMRdZ/bXFixdDLpfjySefFDuKU7EYkWgef/xxHDx4ENu2bRM7ikOUlJRg1qxZ2LRpE1QqldhxnMJisSA1NRUvv/wyAGDo0KE4ePAgVqxY4ZHF6L///S8++ugjrFmzBgMHDkRhYSFmz56NqKgoj3y/1MJkMuG+++6DIAhYvny52HEcIj8/H2+88QYKCgogkUjEjuNUPJQmotDQUMhkMpSXl7fZXl5ejsjISJFSOcfMmTOxfv16bN68GT169BA7jkPk5+ejoqICw4YNg1wuh1wux9atW/Hmm29CLpfDbDaLHdHuunfvjgEDBrTZ1r9/f5w9e1akRI71zDPP4Nlnn8X999+PwYMH46GHHsJTTz2FnJwcsaM5XOvvKG/7/dVaioqLi7Fp0yaPHS366aefUFFRgdjYWNvvr+LiYjz99NOIi4sTO55DsRiJSKFQICUlBbm5ubZtFosFubm5yMjIEDGZ4wiCgJkzZ2LdunX44YcfEB8fL3Ykhxk/fjwOHDiAwsJC2y01NRUPPvggCgsLIZPJxI5od6NGjbps+YVjx46hZ8+eIiVyrMbGRkilbX+NymQyWCwWkRI5T3x8PCIjI9v8/tLpdNi5c6fH/v5qLUXHjx/H999/j5CQELEjOcxDDz2E/fv3t/n9FRUVhWeeeQYbN24UO55D8VCayLKzszF9+nSkpqYiLS0NS5cuhV6vR1ZWltjRHOLxxx/HmjVr8OWXXyIwMNA2F0Gj0cDX11fkdPYVGBh42dwpf39/hISEeOycqqeeegojR47Eyy+/jPvuuw+7du3CO++8g3feeUfsaA4xefJkvPTSS4iNjcXAgQOxd+9evP7663j44YfFjmYXDQ0NOHHihO3706dPo7CwEMHBwYiNjcXs2bPx4osvIjExEfHx8Zg3bx6ioqIwZcoU8UJfh47eb/fu3XHPPfegoKAA69evh9lstv3+Cg4OhkKhECt2l13t5/vr4ufj44PIyEj07dvX2VGdS+zT4kgQ3nrrLSE2NlZQKBRCWlqasGPHDrEjOQyAdm/vvvuu2NGcwtNP1xcEQfj666+FQYMGCUqlUujXr5/wzjvviB3JYXQ6nTBr1iwhNjZWUKlUQq9evYS//vWvgsFgEDuaXWzevLndf6/Tp08XBKHllP158+YJERERglKpFMaPHy8UFRWJG/o6dPR+T58+fcXfX5s3bxY7epdc7ef7a95yur5EEDxkiVYiIiKi68Q5RkRERERWLEZEREREVixGRERERFYsRkRERERWLEZEREREVixGRERERFYsRkRERERWLEZEREREVixGRERERFYsRkTkVZ566incfffd7d63f/9+3H333QgJCYFKpcLAgQOxZMkSNDc3OzklEYmFxYiIvMquXbuQmpp62fYff/wRI0aMgK+vL7788kvs27cPc+bMweuvv467774bFotFhLRE5Gy8VhoReQWj0Qh/f/82oz/p6enYsWMHzGYz+vTpg4yMDHz44YdtHnf06FEMGTIEy5cvxyOPPOLs2ETkZCxGROQVLBYL9uzZg/T0dBQWFiIiIgIqlQpBQUHIy8vDyJEjUVhYiKSkpMsee9ddd0Gv1+O7774TITkRORMPpRGRV5BKpSgrK0NISAiSkpIQGRmJoKAgAMDp06cBAImJie0+NjExEcXFxc6KSkQiYjEiIq+xd+/edkeE1Go1AKC6urrdx9XU1Nj2ISLPxmJERF7jSofKMjIy4OPjg6+//vqy+8xmMzZu3IjRo0c7IyIRiYzFiIi8xoEDB5CcnHzZ9pCQEDz55JN48cUXUVZW1ua+f/zjH6iursZTTz3lpJREJCYWIyLyGhaLBUVFRSgrK0NdXZ1te0NDA5588knExcVh3LhxKCgoAAAsWbIEzz33HN566y0oFAqYzWaxohORk/CsNCLyGh9++CHmzJmDsrIy/PnPf8aSJUsAAH//+9/x/PPP2/abPn063nvvPUgkkjaPP336NOLi4pwZmYicjMWIiIiIyIqH0oiIiIisWIyIiIiIrFiMiIiIiKxYjIiIiIisWIyIiIiIrFiMiIiIiKxYjIiIiIisWIyIiIiIrFiMiIiIiKxYjIiIiIisWIyIiIiIrP4/q5Evq/fz4DgAAAAASUVORK5CYII=", "text/plain": [ "
" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "dynamics_2 = tempo.get_dynamics()\n", "plt.plot(*dynamics_2.expectations(0.5*sigma_z, real=True), label=r'$\\alpha=0.3$')\n", "plt.xlabel(r'$t\\,\\Omega$')\n", "plt.ylabel(r'$\\langle\\sigma_z\\rangle$')\n", "plt.legend()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Finally, we note: to validate the accuracy the result **it vital to check the convergence of such a simulation by varying all three computational parameters!** For this we recommend repeating the same simulation with slightly \"better\" parameters (smaller `dt`, larger `tcut`, smaller `epsrel`) and to consider the difference of the result as an estimate of the upper bound of the accuracy of the simulation." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "-------------------------------------------------" ] } ], "metadata": { "interpreter": { "hash": "3306e98808c0871e8a1685f50cc307ae5b4a4a013844b10634a4efe89132c3fe" }, "kernelspec": { "display_name": "oqupyPR", "language": "python", "name": "oqupypr" }, "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.6.15" }, "widgets": { "application/vnd.jupyter.widget-state+json": { "state": {}, "version_major": 1, "version_minor": 0 } } }, "nbformat": 4, "nbformat_minor": 4 }