{ "cells": [ { "cell_type": "markdown", "source": [ "# Simple Forward Simulation of Burgers Equation with phiflow\n", "\n", "This chapter will give an introduction for how to run _forward_, i.e., regular simulations starting with a given initial state and approximating a later state numerically, and introduce the ΦFlow framework. ΦFlow provides a set of differentiable building blocks that directly interface with deep learning frameworks, and hence is a very good basis for the topics of this book. Before going for deeper and more complicated integrations, this notebook (and the next one), will show how regular simulations can be done with ΦFlow. Later on, we'll show that these simulations can be easily coupled with neural networks.\n", "\n", "The main repository for ΦFlow (in the following \"phiflow\") is [https://github.com/tum-pbs/PhiFlow](https://github.com/tum-pbs/PhiFlow), and additional API documentation and examples can be found at [https://tum-pbs.github.io/PhiFlow/](https://tum-pbs.github.io/PhiFlow/).\n", "\n", "For this jupyter notebook (and all following ones), you can find a _\"[run in colab]\"_ link at the end of the first paragraph (alternatively you can use the launch button at the top of the page). This will load the latest version from the PBDL github repo in a colab notebook that you can execute on the spot: \n", "[[run in colab]](https://colab.research.google.com/github/tum-pbs/pbdl-book/blob/main/overview-burgers-forw.ipynb)\n", "\n", "## Model\n", "\n", "As physical model we'll use Burgers equation.\n", "This equation is a very simple, yet non-linear and non-trivial, model equation that can lead to interesting shock formations. Hence, it's a very good starting point for experiments, and it's 1D version (from equation {eq}`model-burgers1d`) is given by:\n", "\n", "$$\n", " \\frac{\\partial u}{\\partial{t}} + u \\nabla u =\n", " \\nu \\nabla\\cdot \\nabla u\n", "$$ \n", "\n" ], "metadata": {} }, { "cell_type": "markdown", "source": [ "## Importing and loading phiflow\n", "\n", "Let's get some preliminaries out of the way: first we'll import the phiflow library, more specifically the `numpy` operators for fluid flow simulations: `phi.flow` (differentiable versions for a DL framework _X_ are loaded via `phi.X.flow` instead).\n", "\n", "**Note:** Below, the first command with a \"!\" prefix will install the [phiflow python package from GitHub](https://github.com/tum-pbs/PhiFlow) via `pip` in your python environment once you uncomment it. We've assumed that phiflow isn't installed, but if you have already done so, just comment out the first line (the same will hold for all following notebooks)." ], "metadata": {} }, { "cell_type": "code", "execution_count": 1, "source": [ "#!pip install --upgrade --quiet phiflow\n", "!pip install --upgrade --quiet git+https://github.com/tum-pbs/PhiFlow@develop\n", "\n", "from phi.flow import *\n", "\n", "from phi import __version__\n", "print(\"Using phiflow version: {}\".format(phi.__version__))" ], "outputs": [ { "output_type": "stream", "name": "stdout", "text": [ "Using phiflow version: 2.0.0rc2\n" ] } ], "metadata": {} }, { "cell_type": "markdown", "source": [ "Next we can define and initialize the necessary constants (denoted by upper-case names):\n", "our simulation domain will have `N=128` cells as discretization points for the 1D velocity $u$ in a periodic domain $\\Omega$ for the interval $[-1,1]$. We'll use 32 time `STEPS` for a time interval of 1, giving us `DT=1/32`. Additionally, we'll use a viscosity `NU` of $\\nu=0.01/\\pi$.\n", "\n", "We'll also define an initial state given by $-\\text{sin}(\\pi x)$ in the numpy array `INITIAL_NUMPY`, which we'll use to initialize the velocity $u$ in the simulation in the next cell. This initialization will produce a nice shock in the center of our domain. \n", "\n", "Phiflow is object-oriented and centered around field data in the form of grids (internally represented by a tensor object). I.e. you assemble your simulation by constructing a number of grids, and updating them over the course of time steps. \n", "\n", "Phiflow internally works with tensors that have named dimensions. This will be especially handy later on for 2D simulations with additional batch and channel dimensions, but for now we'll simply convert the 1D array into a phiflow tensor that has a single spatial dimension `'x'`.\n" ], "metadata": {} }, { "cell_type": "code", "execution_count": 2, "source": [ "N = 128\n", "STEPS = 32\n", "DT = 1./STEPS\n", "NU = 0.01/np.pi\n", "\n", "# initialization of velocities\n", "INITIAL_NUMPY = np.asarray( [-np.sin(np.pi * x) * 1. for x in np.linspace(-1,1,N)] ) # 1D numpy array\n", "\n", "INITIAL = math.tensor(INITIAL_NUMPY, spatial('x') ) # convert to phiflow tensor\n" ], "outputs": [], "metadata": {} }, { "cell_type": "markdown", "source": [ "\n", "Next, we initialize a 1D `velocity` grid from the `INITIAL` numpy array that was converted into a tensor.\n", "The extent of our domain $\\Omega$ is specifiied via the `bounds` parameter $[-1,1]$, and the grid uses periodic boundary conditions (`extrapolation.PERIODIC`). These two properties are the main difference between a tensor and a grid: the latter has boundary conditions and a physical extent.\n", "\n", "Just to illustrate, we'll also print some info about the velocity object: it's a `phi.math` tensor with a size of 128. Note that the actual grid content is contained in the `values` of the grid. Below we're printing five entries by using the `numpy()` function to convert the content of the phiflow tensor into a numpy array. For tensors with more dimensions, we'd need to specify the order here, e.g., `'y,x,vector'` for a 2D velocity field. (If we'd call `numpy('x,vector')` below, this would convert the 1D array into one with an additional dimension for each 1D velocity component.)" ], "metadata": {} }, { "cell_type": "code", "execution_count": 3, "source": [ "velocity = CenteredGrid(INITIAL, extrapolation.PERIODIC, x=N, bounds=Box[-1:1])\n", "#velocity = CenteredGrid(Noise(), extrapolation.PERIODIC, x=N, bounds=Box[-1:1]) # random init\n", "\n", "print(\"Velocity tensor shape: \" + format( velocity.shape )) # == velocity.values.shape\n", "print(\"Velocity tensor type: \" + format( type(velocity.values) ))\n", "print(\"Velocity tensor entries 10 to 14: \" + format( velocity.values.numpy()[10:15] ))" ], "outputs": [ { "output_type": "stream", "name": "stdout", "text": [ "Velocity tensor shape: (xˢ=128)\n", "Velocity tensor type: \n", "Velocity tensor entries 10 to 14: [0.47480196 0.51774486 0.55942075 0.59972764 0.6385669 ]\n" ] } ], "metadata": {} }, { "cell_type": "markdown", "source": [ "## Running the simulation\n", "\n", "Now we're ready to run the simulation itself. To ccompute the diffusion and advection components of our model equation we can simply call the existing `diffusion` and `semi_lagrangian` operators in phiflow: `diffuse.explicit(u,...)` computes an explicit diffusion step via central differences for the term $\\nu \\nabla\\cdot \\nabla u$ of our model. Next, `advect.semi_lagrangian(f,u)` is used for a stable first-order approximation of the transport of an arbitrary field `f` by a velocity `u`. In our model we have $\\partial u / \\partial{t} + u \\nabla f$, hence we use the `semi_lagrangian` function to transport the velocity with itself in the implementation:" ], "metadata": {} }, { "cell_type": "code", "execution_count": 4, "source": [ "velocities = [velocity]\n", "age = 0.\n", "for i in range(STEPS):\n", " v1 = diffuse.explicit(velocities[-1], NU, DT)\n", " v2 = advect.semi_lagrangian(v1, v1, DT)\n", " age += DT\n", " velocities.append(v2)\n", "\n", "print(\"New velocity content at t={}: {}\".format( age, velocities[-1].values.numpy()[0:5] ))" ], "outputs": [ { "output_type": "stream", "name": "stdout", "text": [ "New velocity content at t=1.0: [0.00274862 0.01272991 0.02360343 0.03478042 0.0460869 ]\n" ] } ], "metadata": {} }, { "cell_type": "markdown", "source": [ "Here we're actually collecting all time steps in the list `velocities`. This is not necessary in general (and could consume lots of memory for long-running sims), but useful here to plot the evolution of the velocity states later on.\n", "\n", "The print statements print a few of the velocity entries, and already show that something is happening in our simulation, but it's difficult to get an intuition for the behavior of the PDE just from these numbers. Hence, let's visualize the states over time to show what is happening.\n", "\n", "## Visualization\n", "\n", "We can visualize this 1D case easily in a graph: the following code shows the initial state in blue, and then times $10/32, 20/32, 1$ in green, cyan and purple. " ], "metadata": {} }, { "cell_type": "code", "execution_count": 5, "source": [ "# get \"velocity.values\" from each phiflow state with a channel dimensions, i.e. \"vector\"\n", "vels = [v.values.numpy('x,vector') for v in velocities] # gives a list of 2D arrays \n", "\n", "import pylab\n", "\n", "fig = pylab.figure().gca()\n", "fig.plot(np.linspace(-1,1,len(vels[ 0].flatten())), vels[ 0].flatten(), lw=2, color='blue', label=\"t=0\")\n", "fig.plot(np.linspace(-1,1,len(vels[10].flatten())), vels[10].flatten(), lw=2, color='green', label=\"t=0.3125\")\n", "fig.plot(np.linspace(-1,1,len(vels[20].flatten())), vels[20].flatten(), lw=2, color='cyan', label=\"t=0.625\")\n", "fig.plot(np.linspace(-1,1,len(vels[32].flatten())), vels[32].flatten(), lw=2, color='purple',label=\"t=1\")\n", "pylab.xlabel('x'); pylab.ylabel('u'); pylab.legend()\n" ], "outputs": [ { "output_type": "execute_result", "data": { "text/plain": [ "" ] }, "metadata": {}, "execution_count": 5 }, { "output_type": "display_data", "data": { "text/plain": [ "
" ], "image/png": "iVBORw0KGgoAAAANSUhEUgAAAZAAAAEGCAYAAABLgMOSAAAAOXRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjMuMywgaHR0cHM6Ly9tYXRwbG90bGliLm9yZy/Il7ecAAAACXBIWXMAAAsTAAALEwEAmpwYAABt+klEQVR4nO2dd3yN1//A3ycbMWOL2GrHiJGEtihKa7XUKkqNTmpTihq/GkW1qD1rFTWqLTXq24rE3nsTW4gVWfee3x/nJrnIzh0J5/163Vfu8zznOedzn3vzfJ5zPktIKdFoNBqNJqU42FsAjUaj0WRMtALRaDQaTarQCkSj0Wg0qUIrEI1Go9GkCq1ANBqNRpMqnOwtgC3JnTu3LFq0qL3F0Gg0mgzF/v3770op8zy//5VSIEWLFmXfvn32FkOj0WgyFEKIy/Ht10tYGo1Go0kVWoFoNBqNJlVoBaLRaDSaVPFK2UA0Gs3LSVRUFMHBwYSHh9tblAyNm5sbnp6eODs7J6u9ViAajSbDExwcTNasWSlatChCCHuLkyGRUhISEkJwcDDFihVL1jl2XcISQswXQtwWQhxL4LgQQvwohDgnhDgihKhqdqyzEOKs6dXZdlJrNJr0Rnh4OB4eHlp5pAEhBB4eHimaxdnbBrIQeDuR442BUqZXD+BnACFELmAEUBOoAYwQQuS0qqQajSZdo5VH2knpNbTrEpaU8l8hRNFEmjQHFkuVcz5ICJFDCFEAeBPYIqW8ByCE2IJSRMutLLLmOW7ehMBAOHkSLl2CGzcgPBwiIyFrVvDwgEKFoEIFqFQJypUDB3s/tphx/dF1goKDOBNyhov3L3LryS2M0ghAwawFKZGzBJXzV+aNom/g4uhiZ2k1mvRFereBFAKumm0Hm/YltP8FhBA9ULMXvLy8rCPlK4TRCLt2werV8PvvcOFCys7PnRsaNIDmzdXLzc06ciZEeHQ42y9uZ8PpDWw6t4nLD+KNj3qB7K7ZafZaMwb4DaBivopWllKT0QgNDWXZsmV89tlnyWp/8eJF2rZtS0hICNWqVWPJkiW4uGS8B5T0rkDSjJRyNjAbwMfHR1fPSiX37sG8eTBjhpppxODuDrVqgbc3FCsGnp6QKRM4OcHjxxASAhcvwtGjsH8/XL0Ky5erV44c0K4d9OkDpUpZT3YpJXuu7WHhoYWsOL6C0PDQ2GPZXLNRs1BNKuatSLGcxSjgXgAnByeM0kjww2DO3jvL9ovbOX7nOEuOLGHp0aV0rdyV0fVGk989v/WE1mQoQkNDmTFjRrIVyKBBg+jTpw9t27blk08+Yd68eXz66adWltIKSCnt+gKKAscSODYLaGe2fRooALQDZiXULqFXtWrVpCZl3Lsn5eDBUmbOLCWoV5EiUg4YIGVQkJTR0cnvy2iU8tQpKadMkbJatbj+hJDygw+kPHbMsrJHREfIJYeXSJ/ZPpKRxL68f/aWI/8ZKfdf3y+jDcn7AKfunJJf/PGFdPzWUTISmf/7/HLvtb2WFViTak6cOGHX8du0aSPd3Nykt7e37N+/f6JtjUaj9PDwkFFRUVJKKXft2iUbNmxoCzGTRXzXEtgn47tHx7fTlq8kFMg7wF+AAGoBe0z7cwEXgZym10UgV1JjaQWSfKKipJw8Wcrs2eNu9A0bSrlxY8qURmIcOSJlt25SOjur/h0cpPzsMynv3Elbv0+jnsqfdv8kPSd7xioNj/Eesu+mvvLIzSNp6vvUnVPyjQVvSEYiM4/NLNefWp82YTUWwfymF/N7tfQrMS5evCjLly8vpZTy4cOH0tvbO97X8ePH5Z07d2SJEiViz71y5UrsuemBDKNAUEbvG0AUyo7xMfAJ8InpuACmA+eBo4CP2bldgXOmV5fkjKcVSPLYs0fKypXj/nHq1VOzDWtx9apSHI6OarxcuaT85Rc1Y0kJ0YZoOXf/XFlwUsFYxVF2Wlk5Z/8cGRYZZjF5I6IjZOe1nSUjkQ7fOsg/z/xpsb41qSM9KZCkeJkUiL29sNolcVwCnydwbD4w3xpyvapER8PYsTB6NBgMUKQITJsG775r3XE9PWH6dPjsM+jdG7Ztgw8/VIb6WbMgb97Ez5dSsvn8ZgZuGcjR20cBqJy/Mt+8/g0tyrTAQVjW7cvF0YUFzReQ3z0/4wPG0+G3DuzvsZ9iOZMXfBVDNLAeyELivuyalCHtbOl89OgRderUiffYsmXLKFu2LKGhoURHR+Pk5ERwcDCFCsXrA5T+iU+rvKwvPQNJmKtXpfT3l7E2iX79pHz82PZyGI1Szp0rZbZsSpaCBaX877+E2x++eVi+tfit2BlHkSlF5LIjy6TBaLC6rAajQTZd1lQyElllZpVkz3KeSilnSCmLSfXDdJVSRlpPzFcCe9tA7t69K728vJLdvlWrVnL58uVSSil79uwpp0+fbi3RUkyGWcKy9UsrkPj57z8p8+aNu2Fv22ZviaS8fFnK2rWVTI6Oyh5jvqT1KOKR7Le5X6xRO/t32eXEgInyadRTm8p5/+l9WWJqCclI5Jd/fplo2wdSynFSynzyxR/nE2sL+pJjbwUipZTt2rWT5cuXT9KILqWU58+fl9WrV5clSpSQrVq1kuHh4TaQMHloBaIVSLKZOzfOiF2/ftoN2JYkMlLK/v1l7Bp0jx5q3++nf5deU7xibRBf/PGFvPvkrt3kPHjjoHT41kE6jXKSp++efuH4AynlGCllThn3Y6wipfxVSpnFtP3QduK+lKQHBfKykBIFko5igjW2REr49lvo1g2iolQsxqZNKtAvveDsDBMnwsqV4OoKs5dfx7NvK5oub8qVB1eoWqAqu7vt5qcmP+GR2cNuclbOX5kulbsQbYzm621fx+5/CIxFuRkOA+4DdYBNwH6gNXGBWNG2FFijsRAvfSCh5kUMBvjiC5g5U6UV+fln6NHD3lIlTOvWkpNOy/h27xfcdgvFIdqdMXXHMODNz3FySB8/4W/f/JZlR5ex5uQatl3bS1Ch6kxCKQ1QiuNbVA4e82xDjqa/BtuJqtFYjPTx36exGQYDdOkCS5aoNCLLl0OLFvaWKmHuPLnDp398ypqTa8ANMl1twtNVM1m2rjBdtkD+dBIMXihbIT7zG8AkYyRv534tdkaRkOKIIeYfUCsQTUZEK5BXCIMBOneGpUshSxb44w944w17S5Uw60+tp8fGHtx+cpusLln54e0faJinCw0CBMeOQb16sGNH0m6+1iYKmAv88uYIEA5EA1UiHjHJNWuCiiOGmBmIXsLSZES0DeQVwWiErl2V8nB3V/aO9Ko8nkQ+ofuG7rRY2YLbT27zZtE3OfLpEbpW6Yqnp+Dff1V235MnVWLGe/fsI6cRWAmUAz4DbgkH8oVehsVvUWvrYOqSuPIAvYSlydhoBfIKICX07QuLF6uZx6ZNULu2vaWKn8M3D+Mzx4e5B+fi6ujKD41+YFunbRTNUTS2TZ48sHUrvPYaHDkCDRvCo0e2lXMLUB1oi0qFUBpYBWwKD4WL21h+bBlPo54m2Y9WIJqMjFYgrwBjxsDUqeDiAuvWgb+/vSV6ESkl0/ZMo+bcmpy6e4pyecqxt/teetfqHW8keb58KmK9eHGV5bdVK1WDxNrsBd4CGgIHgIKoVM/HgVZA5fzeVCtQjdDwUNadWpdkf9oL6+UgJhtvcrl48SI1a9akZMmStGnThsh4frx79uyhcuXKVK5cGW9vb9auXRt7rGvXruTNm5cKFSo8c86AAQMoU6YMlSpVomXLloSGhgJw6dIlMmXKFNvfJ598kroP+jzx+fa+rK9XMQ5kwQIZm6hwzRp7SxM/IWEhstnyZrHR5D1/7ymfRCYvtO7cubggyI4dU54/K7mclVK2knE/phxSBQXGJ+WMPTMkI5H1F9VPst/Spv5OWk7UVxJ7x4GkJBeWlFK2bt36mUj0GTNmvNDmyZMnsRl7r1+/LvPkyRO7/b///U/u37//hTE3b94c22bgwIFy4MCBKZZPx4FoANi+Hbp3V++nTYP33rOvPPFx4MYBqs2uxobTG8jhloNVrVcx892ZZHbOnKzzS5RQzgBZsijPshEjLCvffaAvys6xGnADBgIXgEFAfFK2q9gONyc3tl3cxqXQS4n2r5ewXg4GDx7M+fPnqVy5MgMGDEi0rZSS7du306pVKwA6d+7MunXrXmiXOXNmnJzUHDU8PPyZcrOvv/46uXLleuGchg0bxp5Tq1YtgoODU/uRkoX2wnpJOXUK3n9fJUjs1w/SY62aeQfm8fmfnxNhiMCnoA+rWq96xtaRXHx8YNUqlfRx9GhVNrdt27TJFgX8jHLBvYcyhncBRgGeSZybwy0H75V9j2VHl7Hi2AoG1x6cYFu9hGV5xLfWqY0uRyScpXHcuHEcO3aMQ4cO8ejRIypXrhxvu2XLlpE3b15y5MgRe6P39PTk2rVr8bbfvXs3Xbt25fLlyyxZsiT2nOQwf/582rRpE7t98eJFqlSpQrZs2RgzZkyCCR9TglYgLyEPHqhysaGhKsZj/Hh7S/Qs4dHhfPHnF8w7OA+AHlV7MLXxVNycUl/ftnFjmDwZvvpKxbkULw41aqS8HwlsBPoDZ0z76gKTgCop6Kf5a81ZdnQZf5//O1EFomcgLx9Zs2bl0KFDCR6/e/dusvuqWbMmx48f5+TJk3Tu3JnGjRvjlow60GPHjsXJyYkOHToAUKBAAa5cuYKHhwf79++nRYsWHD9+nGzZsiVblvjQCuQlw2iEDh3gzBmoWFEt6zg6Jn2erbh4/yKtVrXiwI0DuDm58fM7P/NR5Y8s0nevXnD8OMyZAy1bwoEDytieXA4D/YBtpu1SwPdAU5J2x32e+sXqIxAEXA0gLCoswSU5rUAsT2IzBVtgjXTuZcuWxd3dnWPHjuHj45No24ULF7Jx40a2bdsWu+zl6uqKq6srANWqVaNEiRKcOXMmyb6SQiuQl4yRI5VNIFcu5XHl7m5vieL46+xfdPitA/fD71M8Z3HWfLCGyvkrW6x/IZSt5+RJ2LlTLWNt2aLqsyfGDeAbVHEZiSpxOQL4FHBJpSwemT2oVrAa+67v49/L//J2yfgrfuglrJeDrFmz8sjkS57UDASgbt26rF69mrZt27Jo0SKaN2/+QpuLFy9SuHBhnJycuHz5MqdOnaJo0aKJ9rtp0yYmTJjA//73PzJnjntouXPnDrly5cLR0ZELFy5w9uxZihcvnuLP+QLxWdZt9ULV0TmNcqUfHM/xKcAh0+sMEGp2zGB2bENyxnvZvbC2bFG1PBwc1Pv0gtFolGP+N0aKkUIyEvnusnflvbB7Vhvv+nUp8+dXnlkDBiTcLlxKOV5K6S7VD8RJSvmVlDLEQnJ8vfVryUhkn019EmzjZxo7kZInmmRgby8sKS2Tzn39+vXym2++kVJKuXjxYlmuXDnp7e0tq1SpIteuXRt7ftu2bWX+/Pmlk5OTLFSokJw7d66UUsoSJUpIT0/P2BK6PXv2lFJKuXr16mf62rBhQ4KyZYh07qjZ+3mgOOpB7zBQLpH2XwLzzbYfp3TMl1mB3LwpZb586hv99lt7SxPHk8gnss2qNpKRSDFSyNH/G22TYk///htXInd9PGXL/5BSlpJxP45mUsoXE7GnjX8u/iMZiawwo0KCbeqYxt9h4bFfNdKDAnlZyChuvDWAc1LKC1LKSGAF8OI8Lo52qBrqmucwGqFTJ7h1C958E4YOtbdEiqsPrlJnQR1WHl+Ju4s769uuZ9jrwyxeYjY+6tSJcx7o2hWuX1fvzwLvAu+Y3pcBNqNKy5a2sAy+nr5kcc7CsdvHuP7oerxtdDJFTUbGngqkEHDVbDvYtO8FhBBFgGLAdrPdbkKIfUKIICFEi4QGEUL0MLXbd+fOHQuInf74/nv4+29Vy2Pp0vRhNA+8Gkj1OdU5cOMAxXMWJ+jjIJq+1tSmMvTpo9KchIRAh09gkITywB9ANmAycAQVVW4NXJ1ceaOoSji29cLWeNvoZIqajExGCSRsC6yWUpo/qBWRUvoA7YEfhBAl4jtRSjlbSukjpfTJkyePLWS1KUFBcTOOhQuhYEG7igPAokOLeHPRm9x6cou6Reuyp9seyuctb3M5HBxgwUJw/wR2zIQJQsV3dEUZ1PoAzlaWoWFxpZ7+Pv93vMe1F5YmI2NPBXINKGy27WnaFx9teW75Skp5zfT3ArCDlLnpvxSEhkK7dipYsG9feOcd+8pjMBro/3d/Plr/EZGGSD6v/jmbP9xst2qB+4HWBeDxz0BBELth8WmYB6TAuzdN1C9eH4D/rvwX73GtQDQZGXsqkL1AKSFEMSGEC0pJbHi+kRCiDMqzMtBsX04hhKvpfW7AHzhhE6nTEZ9/DpcuqUjs776zryyPIh7RfEVzJgVOwsnBiZnvzGRak2k4O1r7Gf9F7gA9UNlyd6GURb3FIH1hbHN48sR2spTNXZZMTpm48uAK95/ej91vjDZyYvUJytZdROvWq/QSliZDYrc4EClltBDiC5QN0xHlYXVcCDEKZfGPUSZtgRUmT4AYygKzhBBGlBIcJ6V8pRTIunWwbBlkzqyqCrqkNmDBAgQ/DObdZe9y+NZhcmXKxW8f/Ba79m9LDKjMuF8Doagfd29UjIdLa/AZDydOqNnarFm2kcnRwZEKeSuw9/peDt86jG8uXw7MO8CeH/cQeimUrCi7TNTTKMhke2Wr0aQFu9pApJR/SpWQtISUcqxp33Az5YGUcqSUcvBz5+2SUlaUUnqb/s6ztez25N49iMnGPG4clCxpP1kO3jhIzbk1OXzrMKVylSLo4yC7KI+9QC1UYadQlGH8KCqSPDuQKZNStK6uMHu2qoliK7zzeZM9NDtBXwcx2XMyf/f9m9BLoeQqmQvppP4Fo6ONthNIY3Gskc4d4MiRI/j6+lK+fHkqVqxIeHg4YWFhvPPOO5QpU4by5cszeHDc7XHhwoXkyZMnNm373Llz0/zZEiOjGNE1ZvTqpVx269RRy1j2YuOZjdRZUIfrj65Tx6sOgR8HUsqjlE1luI+KGK8J7EMZ0lYDm1AuuuZUqgSjRqn3PXvapghVcFAwxX8qTu+pvYlYHkHko0iKvFGEtuvb8sXpLzBmUbMOrUAyNilVIIMGDaJPnz6cO3eOnDlzMm/ei8/A0dHRfPjhh8ycOZPjx4+zY8cOnJ3V76V///6cOnWKgwcPEhAQwF9//RV7Xps2bTh06BCHDh2iW7duaf9wiaAVSAZj/XrlqpspE8yfrzyN7MG0PdNovqI5T6Ke0KFiB7Z03GJTY7kRWAi8BsxErYEOAE4C75Nw7qq+faFaNbhyBYYMsZJsJvvGPL95zPOdR+TWSKSQBNcIpvu+7ny04yNea/YawkEgnZUZ3RClFUhGxhrp3P/++28qVaqEt7c3AB4eHjg6OpI5c2bq1q0LgIuLC1WrVrV62vaE0LmwMhDpYenKYDTQ7+9+TN09FYARb4xgxBsjnqlVYG2OoJaqAkzbbwDTUbaEpHBygnnzlOPB9OnQpo2ayVmCiEcRHJx3kN1TdxN6KRQAtxxuVOxWkbbGtkTmjOTnKj8/J5BewrI01volJpai0Rrp3M+cOYMQgkaNGnHnzh3atm3LwIEDn2kTGhrK77//Tu/evWP3rVmzhn///ZfSpUszZcoUChcu/HzXFkMrkAxE795w86a64X3xhe3HD4sKo/2a9qw/vR5nB2fmNZtHR++ONhv/ITAS+BFlMM+HsnF0IGU3DW9vGDxYlfrt1g0OH4ZkZMhOWK5rD9n94272z9pPxIMIAHKWyEmtPrWo3LkyLu4u5J6am4uhFzl99/SzMTFagbx0WCqde3R0NDt37mTv3r1kzpyZ+vXrU61aNerXrx97vF27dvTq1Ss2MWLTpk1p164drq6uzJo1i86dO7N9+/bEhkkTWoFkEP74A375RS1dzZtn+6WrkLAQmi5vSmBwIDndcrK2zVqbGcslsBJVGfAGat31C2A0kCOVfQ4bBmvWqMy9336bOjfom4dvEjgpkGPLj2E0KQCv2l749vOldNPSODjGfUne+b25GHqRQzcPPatAnFUbQ5SOBLEU9k3mbrl07p6enrz++uvkzp0bgCZNmnDgwIFYBdKjRw9KlSrFV199FXuOh0fcMnK3bt1emLFYGq1AMgBhYXEzjjFjoJRt7dRcCr3E27+8zemQ03hl92JTh02UzVPWJmOfBj4nrkZHDVSlwKpp7NfVVSlif3+YOBFat4aqyehUSsmFLRfY9f0uLmy5AIBwEJRrXQ7ffr541oy/XmHlfJVZd2odh28dpgMd4g6YZiAGPQPJ0FgjnXujRo2YMGECYWFhuLi48L///Y8+ffoAMGzYMB48ePCCl9WNGzcoUKAAABs2bKBsWev+n2ojegZgzBgVMFi5svLAsiUHbxzEd54vp0NOUylfJQI/DrSJ8ggDhgIVUcojFyrGI5C0K48YfH3V9TQY1FKWIZFJgCHSwKFFh5jpPZNfGv3ChS0XcM7sTI0va/Dl2S9p/WvrBJUHqBkIwOFbh589EKNAtBE9Q+Ph4YG/vz8VKlRI0ogOMH78eCZPnkzJkiUJCQnh448/BtRNf/jw4QDkzJmTvn37Ur16dSpXrkzVqlV55513CA4OZuzYsZw4cYKqVas+4677448/Ur58eby9vfnxxx9ZuHCh1T4zgHg2Pu/lxsfHR+7bt8/eYqSIkyfVmn10NOzaBbVq2W7srRe28t7K93gU+Yi6Reuyts1asrtlt/q4m1BG8oum7Y+BcUBuK4z15ImqoX7liipG9bxbdHhoOPtm7WPPj3t4dF09Ybrnd6dGrxr49PQhU65MyRrn4v2LFP+xOPmy5ONm/5ux+7/xnonTkVtkOdiT/pXzW+xzvWqcPHnS6k/brwrxXUshxH5T7sFn0EtY6Rgp4bPPICoKevSwrfJYemQpH63/iGhjNG0rtGVh84W4OrladcybwFcoewdAJdRylZ8Vx8ySBaZOVSVwhw5VS1l580LopVCCfgji4LyDRD5WQV55yufBr78fFdpVwMk1Zf86RXMUJZtrNm49ucXNxzfJ766UhdBLWJoMjFYg6ZhffoEdO1SadlvlupJS8v2u7xm4VRnf+tbqy8SGE61aw8OIWp4aDDwAMgHfopSJLZJ7NG8OjRvDX3/BNx9fo0GWQE6sOoE0qtl5sfrF8OvvR4lGJVLtriyEoGLeigRcDeDEnROxCiTGiG7URnRNBkQrkHTK/fvQr596P3GiqnFubYzSSN/NfWNjPCY1nERf375WHfMo0JO4TJlNUDEdRa066nNIycDmZyiwKZCCGy9zHHBwcqBi+4r49vMlv4WWlorlLEbA1QAuh16O3adnIJqMjFYg6ZSvv4Y7d1TMR+fO1h8vPDqcTms7serEKlwcXVjUYhFtK7S12nhhwChgEqqYUgFUfEdiUeSWJjo8msNLDhM4KZCQ0yF4AeG4ciVPNX7aU5NcRbNZdDyvbF4AXHlwJXaf0EZ0TQZGK5B0yP79KluskxP8/DNYO8j7UcQjWqxswfaL28nmmo11bdZRt1hdq433F8pIfgmlLD4HxqKSHtqCsLth7J2xlz3T9hB2JwyAbIWzUfWzWnT/uSpnr7hSc4PlPd6K5CgCwOUHcTMQTKlMjHoGosmAaAWSzpBSRZzH/C1v5UJ+d57cocmyJuy7vo/87vnZ1GFTrMuppbmBsmv8atquhLJ91LTKaC8ScjaEoClBHFp4iOinqgJHgaoF8O3nS7nW5XB0duT7csom8s038MEHkN+CjlFFsr+oQPQSliYjo+NA0hmrVkFAAOTJo25i1uTKgyvUWVCHfdf3UTxncXZ22WkV5WFEeVOVQSmPzMBEVPZcWyiPq4FXWdlyJdNem8a+n/cR/TSaUk1K0Wl7J7rv607F9hVxNM0EmjZVlR0fPlTpTiyJV/Z4lrC0Ef2lIKXZeKdNm0bJkiURQiQ7tUl6xK4KRAjxthDitBDinBDihX9XIcRHQog7QohDplc3s2OdhRBnTS8bWAmsz9OnEBODNGYMZLfims6pu6fwn+/P6ZDTVMxbkZ1ddlIiV7xl5dPEEVS5yM9QuayaAMeB/ljXw0oaJWc2nmFBnQXM95vPqXWncHR2pMrHVfjs+Ge0/6M9xeoWe8GrSgjl1uviAosWqeVES2GuQIxSzThiZiB6CStjk1IF4u/vz9atWylSpIgVpbI+dlvCEkI4ohxuGgDBwF4hxIZ4KguulFJ+8dy5uYARgA8q9c1+07n3ycBMmqQC2ipVAlNgqlXYe20vjZc2JuRpCP6F/fm93e/kzJTTomM8Ic5IbsB2RnJDpIGjy46ya+Iu7py4A6iMuD6f+VDzy5q453dPso8SJZT94/vvVfr3HTssY4fK4pKF3JlzczfsLrce36JA1gI4ODlgRCuQjI55OvcGDRowceLERNtXqVLFRpJZF3vaQGoA56SUFwCEECuA5iSvtnkjYIuU8p7p3C3A28ByK8lqda5fj4v1+OEHcHS0zjjbLmyjxcoWPI58TOOSjVn9wWoyO2e26BibgU+wrZE84mEE++fsJ2hKEI+umXISFcqKb19fqnavimvWlAVBDh0KCxfCv//C2rXw3nuWkbNI9iLcDbvL5QeXKZC1ACLGiK69sCzGt+Jbq/Q7Qo5I8FhK0rmXK1fOKvLZA3sqkELAVbPtYOJfEn9fCPE6cAboI6W8msC5L6azzEAMGaKSJrZsCXWt5AD128nfaLemHZGGSNpXbM/C5gtxdrTcQlII0AdYYtq2hZH88c3H7P5xN3tn7I1NpZ6nXB78BvpRsV1FHF1Sp4lz5FBZej//HAYOVHYRVwsE4ntl92L/jf1ceXCFWp61cDAtYUk9A3lpSE4yxZeF9O6F9TuwXEoZIYToCSwC6qWkAyFED6AHgJeXl+UltAB79sDixWrdPYmZb6qZd2AePTb2wCiNfFH9C6Y2nmqx6HIJrAB6A3cAN1Tdjr5Yz84RciaEXd/v4vCiwxgilQHaq44X/gP9KdWkFMIh7WtOPXqo/FgnT6q/MYGdaSHWE8sUTCiclJzaiG45Epsp2IKk0rnrGYhluAaYl8ryNO2LRUoZYrY5F5hgdu6bz527I75BpJSzUQ/C+Pj4pLvMkVKCKUMzX32l1t8tzcSAibGpSb5981u+ef0bi1UQvIKqSf6naftN1MW2Vsb5a3uuETA+gJNrTyrNJaBMizL4DfSjsK9lK685OSm7VJMmMHo0dOqkvOPSwvOxIA46DuSlIKXp3F8W7OmFtRcoJYQoJoRwAdoCG8wbCCEKmG02Q5W8BrXM3lAIkVMIkRNoaNqX4fjtN5VlN29ete5uSaSUDNs+LFZ5/NT4J4a/MdwiysMATEOVkf0TZd+YC2zH8spDSsnZv86y8M2FzK05l5O/nYz1qPr8xOe0WdvG4sojhsaNoVEjePBALWmllRhPrFgFopewXgpSms79xx9/xNPTk+DgYCpVqkS3bt2SPCc9YrcZiJQyWgjxBerG7wjMl1IeF0KMAvZJKTcAvYQQzVDZLu4BH5nOvSeEGI1SQgCjYgzqGYmoKGX7ABgxArJZMHOGlJK+m/vyw+4fcBSOLGyxkA8rfWiRvo8D3YnLX/U+8BPK08qSGKIMHFtxjF0Td3H76G0AXLO54vOpDzV71SRrwawWHjF+Jk2CLVtg5kyVHTktKxAxS1gxsSAOpjgQqY3oGZ5ly5Ylu22vXr3oZeviPlbArjYQKeWfxK1+xOwbbvZ+CDAkgXPnA/OtKqCVmTsXzp5VFQa7d7dcvwajgU//+JQ5B+bg7ODMylYraVm2ZZr7jQC+A/4PiEIpjOlA2nt+lsjHkRyYe4DAyYE8vPoQAPcC7tTqU4tqParhlj0NBcxTQfnyyh4yc6ZS+OvXp76v2CWsUD0D0WR80rsR/aXl8eO4JZH/+z9wtpC1OdoYTed1nVl2dBluTm6sbbOWt0u+neZ+A4FuxPlY90QVecqR5p7jeHL7Cbt/2s3e6XsJvx8OQO4yufEb4EfFDhVTXIPDkowYAUuWwIYNKlOAv3/q+vHI5EFm58w8iHjAg/AHcQpEG9E1GRCtQOzEpElw6xbUrAnvv2+ZPiOiI2i3ph1rT63F3cWdje028kbRN9LU5yPga9RMQ6LsG3OAtPX6LPfO3yNwUiCHFhwiOlzlqPL09cR/kD+vNX3NIh5VaSV/fuWFNWoUDBoE//2XuuBCIQRe2b04dfcUlx9cjk2homcgaUdKaTHnkFeVlFao1QrEDty6FeeuO3GiZaKcw6LCeG/le2w+v5kcbjnY1GETNT3TFoHxB8rD6irqhzIQ+AblpmsJru+/rjyq1pyMLd5Uumlp/Af641U7/blc9+sHM2aoGcjGjSpvVmookr0Ip+6e4sqDKzg4KTuOViBpw83NjZCQEDw8PLQSSSVSSkJCQnBzS/5/uFYgdmDUKFWLu2lTVe8jrTyKeMS7y9/l38v/kidzHrZ03JKmpIi3UTEdK0zbPigPK0ukWZRScmHLBQLGB3Bxu6p67uDsgHcnb/wG+JGnXBr9ZK1ItmwqwWXv3soW0qRJ6jIGmMeCODtXArQRPa3EeDTduXPH3qJkaNzc3PD09Ex2e61AbMzZszB7Njg4wLhxae/v/tP7vL30bfZc20PBrAXZ2nErZfOUTVVfEhVF3gfl8pYJGAP0Iu0/FGO0keO/HidgQgC3Dt8CwMXdhWo9q1Hrq1pk87Rs8SZr0bMnTJkCx48rm8hHH6W8D/NYkNecKqudegaSJpydnSlWrJi9xXjl0ArExgwdCtHRKlliWgNSbz+5TcMlDTl86zBFcxRlW6dtFM9ZPFV9XUQZxreYthsAM4HU9RZH1NMoDi04xK6Juwi9FApAlnxZqNm7JtU/rY5bDtt6VKUVV1eVKfnDD9VspE0byJQpZX0UzqZiVq4+vEo57YWlycBoBWJDDh5U9T7c3NIelHbt4TXeWvIWp+6eorRHabZ12oZntuRPPWMwAFNRto0wICcwBehE2rLmhoeGs/fnvez+YTdPbj8BIFfJXPgN9MO7ozdObhn3p9eunbJdHT4M06dD//4pO79g1oIA3Hx8M9aIjvbC0mRAMu5/cQYkpkDUZ59BoTSkfrx4/yL1F9fnYuhFKuatyJaOW8jnni/F/RwDuhIXjdkW+AFIeU9xPLrxiKAfgtj38z4iH0UCquqf/2B/yr5XFgfHjF/DLGb5sXFj5YLdrZtKvphc8rurMoc3Ht3QcSCaDI1WIDYiMBD++AOyZElbpbvTd09Tf3F9rj26RvWC1dn04SZyZcqVoj4iUQGBY1EBgZ6oioHvpl4s7p27x67vd3Fo4SEMEeppuli9YvgP9qf4W8VfOs+YRo3gzTdVrZDJk5VjRHIpkFXF7N94fAOn3CaFqo3omgyIViA2Ytgw9ferr1KfkO/IrSM0WNKA209uU8erDhvbbySba8qMz/tQs46jpu1PgPFAak3YNw7eIGB8ACdWnYh1xS3Tsgz+g/zxrJnyJbWMghDKFlK7tjKq9+oFuXMn79ycbjlxcXThYcRDjA6mpSs9A9FkQLQCsQHbt6tXjhwpXy+PYe+1vTT6pRH3w+/TsERD1rZZm6JCUE9RKda/R9UoL4FyzX0zFbJIKbn872UCxgVwbtM5QKXk8O5scsUtm35dcS2Jvz+8/TZs2qRsIuPHJ+88IQT53fNz5cEVwoxhaqdWIJoMiFYgVkbKONtH//4pWyuPIfBqIG8vfZuHEQ9p/lpzVrZaiatT8qsb/Qd8DJxFpV/uC4wGUlqHMKbO+M5xOwkODAbAObMzVXtUxbevL9kLW7PmYPpk1CilQKZNU+Vv8yXTgBSjQJ4alYOBXsLSZES0ArEyf/2l0rXnzq2WOVLKzis7aby0MY8jH9O6XGuWvrc02VUEH6EyUU43bZdDZZ9MaXy6IcrAseXHCBgfEFtnPFOuTNToVYMaX9Qgs4dlS+JmJKpXh2bNVI6scePUclZyKOCu7CCPo1UNCaFnIJoMiFYgVkTKONvHkCGQNYXZx3dc2sE7y94hLCqM9hXbs6jFIpwckveV/Y1KuX4F9SUPAYYCKanKGhUWxYF5Bwj8PpAHVx4AkM0zG779fKnarSou7i4p+jwvK6NGKQXy888q3UlyAnnjFMhDVbVRu/FqMiBagViR335TsR8FCsCnn6bs3K0XttJseTOeRj+lk3cn5jebj6ND0nkz7qOWqBaatquiZh0pSUPy9P5T9kzbw54f9xB2V63R5y6TG/9B/lRsn/o64y8r3t7QurWK8fm//1P5spIixpX3keEBudAzEE3GRCsQK2EwwHBTZZNhw1IWrbzp3CZarGhBhCGCblW6MavprGTVL18LfAbcRM00vgX6kfwv+eG1hwRNCWL/rP1EPlYxHAWrF6T2kNqUaV4mXWTFTa+MHAmrV6saL4MGQZEiibePceV9EHmfXLhoI7omQ2JXBSKEeBsVCO0IzJVSjnvueF9UGYpo4A7QVUp52XTMQJw36hUpZTObCZ4MVqyAEyfUjSQl1Sp/P/07rVa1ItIQyac+nzKtybQklcdt4EvgV9O2PzAPeC2ZY4acCSFgQgCHFx/GaDLmFm9QnNpDalP0zaIvXQyHNShXDtq3h6VLVf30uXMTbx8zAwmNvg/kQ2gjuiYDYjcFIoRwRNl3GwDBwF4hxAYp5QmzZgcBHyllmBDiU2AC0MZ07KmUsrItZU4uUVGqABGovy7JNBWsPbmWNqvbEGWMoleNXvzw9g+J3rwlsAyVOTcEyIIKEPyc5BW7v77/OgHjAjix5oTqTEC51uXwH+RPwWoFkye0JpYRI9SDw8KFKli0ZMmE28bYQO5HhAD59BKWJkNizxlIDeCclPICgBBiBdCcuKJ3SCn/MWsfBFimqLeVWbQIzp+H0qWhY8fknfPr8V9pv6Y9Bmmgn28/JjaYmKjyCEYFAf5h2n4LVeipaBLjSCm59M8ldn63kwtbLwCmdOqdvfEf4I9HaY/kCax5gVKloHNnmD9f5TpbsiThtjFLWPci7gLaBqLJmNhTgRRC1SqKIZjEPUw/Bv4y23YTQuxDLW+Nk1Kui+8kIUQPoAeAl5f1ixRFRMSltfj2W3BKxhVednQZHdd2xCiNDKk9hLH1xiaoPCRKUQwAHgLZgclAFxJPfiiNklPrTrFz3E6u770OmNKpf2JKp14oY6RTT+98841SHEuXKs+7hDIu582SF4B7kSYFor2wNBmQDGFEF0J8iKprZF5JtYiU8poQojiwXQhxVEp5/vlzpZSzgdkAPj4+KavXmAoWLoSrV6F8efjgg6TbLzq0iK4bumKURoa/PpyRb45MUHmcR7nmxkzLmgMzgMQWmwyRBo4sPULA+ABCTocAkDl3ZpVO/fPqZMqZwlzkmkQpWlSl6p85UxnWf/01/nYuji7kzpwbQ6hyVtAzEE1GxJ4K5BpQ2Gzb07TvGYQQb6FCGN6QUkbE7JdSXjP9vSCE2AFUQd1j7UZkpHLjBOWB5ZCEIWLugbn0+L0HEsnouqMZ9vqweNsZgB9RF+EpkBuYBnxAwrOOyMeR7J+zn8BJgTy6poLVsntlx2+AH1W6VsE5c/KCETUpZ+hQWLBAufUePqzcfOMjv3t+bomnADhoI7omA2JPBbIXKCWEKIZSHG2B9uYNhBBVgFnA21LK22b7cwJhUsoIIURulOPRBJtJngCLFsGVK1C2LLz/fuJtf977M5/9+RkA498az0D/gfG2O4lKfhhk2m6HcltLKNtU2N0wFcPx0x6e3lM3pzzl8uA/2J8KbSvE1Z/QWA1PT/jkE5g6VS1nrlkTf7sC7gW4IU4DegaiyZjYTYFIKaOFEF8Am1FuvPOllMeFEKOAfVLKDcBEwB1YZVrWiXHXLQvMEkIYUQ5H457z3rI5UVFxs49vvkm8VvaPu3+k96beAExuOJk+vn1e7A+lEUeh0q8XRFUIbJpAnw+uPCBwciAH5hwgKiwKAM9antQeUpvS75bWMRw2ZtAgtYz1229w9ChUrPhimwJZC2AUxwGtQDQZE7vaQKSUfwJ/PrdvuNn7txI4bxcQz7+k/ViyBC5dgjJlErd9TNo1if5bVErenxr/xBc1vnihzUHUrOOQabsbSpPmiKe/OyfvsGvCLo78cgSj6SZUsnFJag+ujVcdLx3DYScKFIAePeCnn1RcSHy2kPxZ8mMQygbioI3omgxIhjCip3eiomDsWPV+2LCEZx/f/fcdX2//GoCZ78ykp0/PZ46Ho2YcE1B2j2Ioj6v68fR1bc81dn63k1PrTgEgHAQV2lbAf5A/+SvnT/uH0qSZQYNg1iwVoX78uHKsMEfNQKIBcNAzEE0GRCsQC7BsGVy4oOIA2rSJv82o/41ixI4RCARzm82la5WuzxzfhfJTPoUyjPdGVQzMYtZGSsmFLRfYOW4nl/65BICjqyOVu1TGr78fuUqkrDKhxroUKqSyEMyYoYpPLV/+7PH87vkxopYbHaKMMbGcGk2GQSuQNBIdrW4OoGYfz8d9SCkZ/s9wxvw3BgfhwMLmC+noHRdd+ATlXfUjKsbjNVTyQz+zPowGIyd/O0nAuABuHLgBgEtWF6p/Vp2avWuStUAK0/xqbMbgwTBnDqxcqTzzypaNO1bAvUBsRUKHaCNGlDFQo8koaAWSRpYvh3PnVNqK9u2fPSalZPDWwUzYNQFH4ciSlktoV7Fd7PHtKPvGRdSNYyAwHHAzHY+OiObw4sPsmrCLe+fuAZAlbxZq9amFzyc+uOVwQ5O+KVw4Li5kzBgVYBhDgawFMDqopSuHaCMGtALRZCy0AkkDBkPc7GPo0GdnH1JK+v3djylBU3BycGL5+8tpVa4VAA9QymK2qa03atZR1bQd8SiC/bP2Ezg5kMc3HgOQo1gO/Ab4Ufmjyjhn0jEcGYnBg2HePJUna/hweM2U5TK/e34MjmoG4hhlQJvRNRkNrUDSwMqVcOYMFCsGHTrE7TdKI73+6sX0vdNxdnBmVetVNC/THFC5WHqg8rY4o2Ycg0zvn9x5wu6pu9k7fS/hoeEA5KuUD//B/pRvXR4Hp+SkSNSkN4oUgY8+UktZY8fC4sVqf1aXrLi6qBJfDtFGou0nokaTKrQCSSUGg3LPBDX7cDZNCozSyKcbP2X2gdm4OLrw2we/8U7pd7gH9AFM9w6qo2YdFYDQS6HsmrSLg/MOEv1U3Ua86nhRe3BtSjYuqV1xXwK+/lpFpy9dquKESpUCIQR5s6qcWEJCtFGCjtfRZCC0Akklq1bBqVMq91GnTmqfwWig++/dWXBoAW5Obqxrs45GJRuxFvgUuIWyb4wGvgLuHbvN2vEBHF1+FGlQabpKv1sa/8H+ePlbP/GjxnbE/E7mz1cBpwsWqP153fNicHbAMcpIVJQBXPW/pCbjoH+tqcBojJt9fP21mn0YjAa6rO/CkiNLyOSUid/b/U7F4vVpQ1yhpzrAXCDTrqusGreTM7+fAUA4Cip9WAm/gX7kq5jPDp9IYwuGDlXpbpYsUR57JUqorLxGJ5MCiTamrGi9RmNntAJJBWvWqGqDXl6q/kO0MZqOazuy4tgKsjhnYWP7P7he9A3aEFfoaZyUNPjrHAHjdnLlvysAOLk5UeXjKvj19yNH0Rx2/EQaW1C8uKoPs3AhfPedqlqYN3NejE5q2SpSJ1TUZDC0AkkhRmNcvY8hQ0A4RtF2dTvWnFxDVpesLO60jcmFqvO7qf1b0UaG/HqcM+MDWHHkFgCu2V2p8UUNavaqSZa8WeIfSPNS8vXXyoi+aJGaheRzz4fB5BwRraPRNRkMrUBSyNq1cOyYyrjavmMErVe1Yf3p9WRzzc6XPQ7wUa7iPAByPY1i+IJDMHEX/10KBcC9gDu+fX2p1qMartn0WsWrSKlSymNvyRJlCynXNS/XnFSVAq1ANBkNrUBSgPnso//gcNqvb8UfZ/8gW96KlO3yH2PdsuN2/yndZ+ylxNTdhN4JAyBXqVz4D/SnUsdKOGkj6SvP0KHKG2vhQvi+XV6MTqoMTqROqKjJYOi7WQrYsAGOHIH8hcPY6N6SrWe3kMW3P9FvjePEzTCaTvkbn5n7kY8jCQcKVCtA7SG1KdOiDA6OOoZDo3jtNWjbVuVQ27IuLyWdVIlhPQPRZDS0AkkmUppmH64PydKzKVsfXMO5607cnr6GX88/qLL4CA6RBiRQ/K3i+A/2p1i9YjqGQxMvw4apNDib1uTlM2f1G4nWRnRNBkMrkGTy++9w8NQ9nHo25bzvmxTM/DH+ffdTbs0WhCmNarnW5fAf5E/BaolVKddoVFLFDz6AlX/kxVBAKZAoPQPRZDDsqkCEEG+jKrQ6AnOllOOeO+6KCt6uhvKIbSOlvGQ6NgSVAd0A9JJSbraWnFLCsHG3cBj9Na9lG0HNb4Lx2rUEAAcXRyp39savvx8epT2sJYLmJeSbb2DlqtwYTRkUIyIj7SuQRpNC7KZAhBCOwHSgASo11F4hxIbnStN+DNyXUpYUQrQFxgNthBDlUDXUy6OqvW4VQpSWUlrFCjl1xXlyeG/lqynlyXYtEADH7K7U6F4V3z6+ZC2o06lrUk758tD6fSeMJ5R9LPTJQ6CQfYXSvHSEh8PDsEjy5nKxeN/JUiBCiOHx7ZdSjkrD2DWAc1LKC6YxVgDNAXMF0hwYaXq/GpgmlFGhObBCShkBXBRCnDP1F5gGeeJl++mzhHRZSf0IpZuMZXLRuLcvVT+shIu75b8QzavFsGEwraNSIFduPbCzNJqXjejwaHoM/hPngCt4tKnNhP6VLdp/cmcgT8zeuwHvAifTOHYh4KrZdjBQM6E2UspoIcQDwMO0P+i5c+N9dBNC9EAlwMXLK+X5pTxdCnPVOweRmQRtv3qdps0raMO4xmJUqkRsJPqOnQ/p856dBdK8FDy++ZhNU4M4OHs/Je6pzN7HKxwCKlt0nGQpECnlJPNtIcT3gNVsDpZESjkbU+kNHx8fmdLzSxdzY/gfnbh2UVCnul6q0lge4aKMIMdOP+TWLcin06FpUsm9c/cI+H4XBxYegggDLsCNynk5XM2ZzZPbJXV6ikmtDSQz4JnGsa8Bhc22PU374msTLIRwArKjjOnJOddiFM+djeK5rdW75lVHOJtihBwf8/33MHGifeXRZDxuHLxBwPgAjq86AUb1nHyyZRkCWjwk+I9rrGw4gmzZLL9yklwbyFFUyW5QHlN5gLTYPwD2AqWEEMVQN/+2wHNFYdkAdEbZNloB26WUUgixAVgmhJiMMqKXAvakUR6Nxi4IUy4sB5cnzJgBAwdCnjx2FkqT7pFScmnHJQLGBXD+7/MAGJwdONLZm4BPinP37EBYUJky1/+PVsuts+ye3BnIu2bvo4FbUso0FVAz2TS+QC2FOQLzpZTHhRCjgH1Syg3APGCJyUh+D6VkMLX7FWVwjwY+t5YHlkZjdVzUv2GO3E84FwaTJsG4cUmco3llkUbJqfWnCBgXwLU9auElOosze3pWI6iPLw9vLoctHcmxtw+h20byzVKBg5USYSTXBnLZGoNLKf8E/nxu33Cz9+FA6wTOHQuMtYZcGo0tEU7KBpLT4ykA06ZB//6QWy+baswwRBo4svQIAeMDCDkdAoD0yMT/etdk9+c1yJIliocrWsClf2iZfTRr1w6jVClo08Z6MulIdI3GzjiYZiBGwnn7bdi0CX74AcaMsa9cmvRB5ONI9s/eT+DkQB5dewSAm1d2dvb3ZXPXKkRncaHuzUNsn+kPUWGMqzeBOV0GAMpN3NHRerJpBaLR2BkHZ/Uf/jTiKaOHKwXy44/Qty/kymVn4TR248mdJ+z5aQ97pu0h/L5yxfUon4fLg/wZ3bYCBmdHygINj61k6pq2APzQ6AdynunN+fNQsiS0f96qbGG0AtFo7IyjszMAERER+PpCgwawZQtMnQrffmtn4TQ2J/RyKIGTAjkw9wDRT5WpubBfYbIOqc2gJqW46CBwBIYCmXZNYtiW/gBMbzKdHlU+o1wn1c+wYeBk5Tu8zjGu0dgZB2fTEpZB8iTyCcNNVsCpUyE01H5yaWzL7WO3WdtxLT+W+JE9P+0h+mk0pd4pRev/uhAU0JUP3i3NRQdBZZQLq9u/Yxi2pT8CwZymc/is+mesWAFnz0KJEqpwmbXRMxCNxt6Y4kAcjA7cfnKb2rWLUa8ebN+ulrKGx5tISPOycCXgCgHjAjiz8QwAwlFQsUNF/Af6s69SPhqi4hxcgOHAACkZs2MEo/8djYNwYH6z+XSu3BmDIc5uNnSo9WcfoBWIRmN3YuJAHI2O3H5ym2I5izF8uFIgU6ZA796QPbudhdRYFCklZ/88S8C4AK7svAKAk5sTVbpVwa+fH4aiOegNLDW1rwnMB8pKyZBtQxgfMB5H4cjilotpX1EZOlauhNOnoVgx+PBD23wOrUA0GjsjYmcgjtx6cguAN96A11+Hf/9Vbr1Dh9pTQo2lMEYbObbyGAHjA7h99DYAbjncqP5FdWp+WZMsebOwGvgcuA1kAsYAvQEHKen3dz+mBE3BycGJZe8to3V5FeVgMMDo0WqMoUPBZFazOlqBaDR2JjYS3ejArce3YvePGAH168PkydCrF2TVqdgyLFFhURycf5Bd3+/iwWWVdTlrwazU6luLaj2q4ZrVlZtAJ+A30zlvAnOAkoBRGvnyr15M3zsdZwdnVrVeRfMyzWP7X7UKTp2CokWhUyfbfS6tQDQaOyNMbrwOBofYGQhA3brg7w8BATB9OgwebC8JNanl6f2n7J2+l91TdxN2NwwAj9Ie+A30o9KHlXBydUKiquZ9BdwHsgITUCnEHVDK45ONnzDnwBxcHV35rc1vNCnVJHYMozFu9vH117abfYBWIBqN3TGfgdx8fDNuv1CzkIYNVXqTL74Ad3d7SalJCQ+vPSRoShD7Z+0n8rGqNFnQpyC1h9Tmteav4eBoqgED9AQ2mc57G5gFxBSeMBgNdPu9GwsPLcTNyY31bdfTsETDZ8ZavRpOnAAvL+jc2QYfzgytQDQaO+MQY0Q3OHDj8Y1njr31FtSqBUFB8PPPMGCAPSTUJJe7p++ya+IuDi8+jDFK1bgv/lZxag+pTdG6RWNrCRlRy1MDgEdATuAHoCMQk/Yw2hhN53WdWXZ0GZmdM/N7u9+pV6zeM+MZjTDKlNb266/BxcY17rQC0WjsjDBz473x6FkFIoRy423SRKV5//xzyJzZHlJqEuP6vuvsHLeTk7+dVHnLBZRrXQ7/Qf4UrFbwmbbngW7ADtN2S2AGkN+sTZQhig6/dWDViVW4u7jzZ/s/qVOkzgvj/vYbHD8OhQtDly7W+GSJoxWIRmNnHGKXsMQzS1gxvP02VK8Oe/fCrFnQp4+tJdTEh5SSi9svEjAugAtbLwDg4OyAd2dv/Af441Ha45n2BuBHVAT5U1RNjOmoOhXmydYjoiNou6Yt606tI5trNjZ12IRvYd8XxjeffQwZYvvZB2gFotHYHXMj+o3HN5BSPlM2OWYW0rQpTJgAn3wCmTLZS1qN0WDk1DqVTv36vusAuLi7UO2TatT6qhbZCmV74ZyTQFfi6nB3QC1ZPZ9wOTw6nFa/tuKPs3+Qwy0Hf3/4N9ULVY9XjrVr4ehR8PSErl0t89lSilYgGo2diZmBOEknwqPDeRjxkOxuz0YOvvMOVK0KBw7AnDnKrVdjWwyRBo78coSACXHp1DPnzkzNr2pS/bPqZMr5olaPQnlUjQIigULATJ4tsBTD06intFjZgr/P/41HJg+2dNxClQJV4pXFfPYxeDC4uqb986UGrUA0GjsT44XlilqDuPH4xgsKJGYW0qKFKjbVvbuehdiKiEcRHJhz4Jl06tmLZMdvgB9VulTBOXP8frMHULOOw6btbsD3qLrcz/Mk8glNlzfln0v/kCdzHrZ12kbFfBUTlGn9ejhyBAoWhI8/TsOHSyN2USBCiFzASqAocAn4QEp5/7k2lYGfgWyo5cOxUsqVpmMLgTeAB6bmH0kpD1lfco3G8jg6OyAxUyCPblAmd5kX2jVrBpUrw6FDyhby1Ve2lPLVI7506nnK56H24NqUb1MeR+f4C22EA98CE1E3rmLAXKBevK3hUcQj3ln2Dv9d+Y/87vnZ1mkb5fKUS1AuKZ+dfbi5pe7zWQJ7zUAGA9uklOOEEINN24OeaxMGdJJSnhVCFAT2CyE2SylDTccHSClX205kjcY6ODg5YABcpHqSjc+QDmoWMmqUUiTffadmIVmy2FDQV4TQy6EETg7kwJxn06nXHlKbUk1KIRwSri++E/gYOIMyjH+FSkWS0Nf0IPwBjZc2JjA4kEJZC7G983ZKe5ROVL4NG9RDRIEC6jdgT+ylQJqjIvUBFqE82p5RIFLKM2bvrwshbqMcF0JtIqFGYyNiFIizUf07Ph8LYs6770KNGrBnj8qRNej5xy5Nqrl9/Da7Juzi6LKjGKNVDEepd0pRe3BtvGp7JXruY2AIyqtKAmVRyQ9rJXLO/af3afRLI/Ze34tXdi+2d9pOiVwlEh3HaFTBpaC+e3vOPsB+CiSflDLmv+QmkC+xxkKIGqhsxufNdo8VQgwHtgGDpZQRCZzbA5UVAC+vxH8EGo09iKlI6GRUfxOagYCahYweDY0aKY+sTz+FbC86/WhSwNXAqwSMC+D0htMACAdBxfYV8R/kT75Kid6aAPgbdYO5jLqhDgaGAYnZte+G3aXhkoYcvHmQYjmKsb3zdormKJrkWGvWwOHDUKgQ9OyZZHOrYzUFIoTYyrOxMTE8k1dUSimFEDKRfgoAS4DOUkqjafcQlOJxAWajZi+j4jtfSjnb1AYfH58Ex9Fo7IVjrBeWUiCJzUBAVSysXRt27lS103W9kJQjpeT85vPsHLeTy/+7DKh06pW7Vsavnx85i+dMso/7QF9goWm7KmrW4Z3EeTcf3+StxW9x/M5xSuYqyfZO2ymcvXCS4xkMcd/1N9/Yf/YBVlQgUsq3EjomhLglhCggpbxhUhC3E2iXDfgDGCqljHGhxmz2EiGEWAD0t6DoGo1NcTBFojsZ1N/no9GfRwhVOOjNN1Wm3i+/hJxJ3+80qHTqJ1afYOe4ndw6rBJXumZ3pfrn1anZqybu+ZKXbGwt8BnqKdYVZTTvR9I31OCHwdRfXJ8zIWcol6ccWztupUDWAskac+lSlXG3WDH7RJ3Hh72WsDYAnYFxpr/rn28ghHBBfU+LnzeWmykfAbQAjlldYo3GSjjGFpRSxtnElrBieOMNlep92zaVaDGmEp0mfqLDozm08BC7Ju7i/gXl8Ome351afWrh84kPrtmSF0hxC/gSWGXaro3ysHotGedeCr1EvUX1uBh6Ee983mzpuIU8WfIka9zISBg5Ur0fOdI+UefxYS8FMg74VQjxMWrp8AMAIYQP8ImUsptp3+uAhxDiI9N5Me66S4UQeVCODoeAT2wqvUZjQWIViEEpkKSWsGIYPVopkB9+UFUL8yTvXvRKEf4gnH0z9xE0JYgnt54AkLNETvwH+uPdyRsnt+TdAiXwC8qr6h7Kq2o88Ckq5XpSnA05S73F9Qh+GEz1gtXZ9OEmcmXKlezPsWABXLwIZcrYptZ5crGLApFShgD149m/DxVvg5TyF9R3Ft/5CblUazQZjhgjujCAo3Dk3tN7RERH4OqU+FOxry80bgx//aUM6hMn2kLajMHjm48JmhrEvhn7iHio/GvyV8lP7cG1Kft+2dh06snhKirl+l+m7YYoo2qRZJ5/4s4J6i+uz83HN/Ev7M8f7f94IVA0McLD4+p9jBoFjvGHn9gFHYmu0dgZJ9MMhGhJPvd8XH90nVtPbuGVPWmvwdGjlQKZPh369lWxAa8y987fY9f3uzi04BCGCAMAResWpfbg2hRvUPyZHGNJYUQpioGolOs5gCmoNffk9nL45mHeWvIWd8PuUrdoXTa024C7S8qKusycCdeugbc3vP9+ik61OlqBaDR2xtFkRCfKQAH3Alx/dJ0bj24kS4FUq6bSm6xbp4ILf/zRqqKmW24evknA+ACOrzyONCpnyzItyuA/yB/PWp4p7u8s0B34n2n7PVSMR3xupQmx99peGv3SiPvh93m75Nv89sFvZHJOWf6Zx4/V9wrqYcEh+RMnm6AViEZjZxxjZyBG8rurW1RyDOkxfPutUiCzZqmCU4WT9gh9KZBScuW/K+wct5Nzf50DVFCmdydv/Ab6kadsyo1C0agsud+gUpLkJS7lekoIuBJA46WNeRT5iOavNWdlq5VJLknGx7RpcPs21KypgkjTG1qBaDR2xlyBFHBXa1DJNaQDVKoEH3wAv/4KY8eqJY+XGWmUnNl4hp3jdhIcGAyAc2Znqvaoim9fX7IXTr59wZyjqDQke03bHVFLVh4JnhE//1z8h6bLm/Ik6gltyrdhScslODumvFB5aKiybYHyskvB6pvN0ApEo7EzTjFJ+aJSNwMB5dq5ejXMmwcDB0Lx4hYWMh1giDJwbMUxAsYHcOf4HQAy5cpEjS9rUOPLGmT2SF2pxkjg/0yvKMATVZe8SSr62nRuEy1XtiQ8OpxO3p2Y32w+jg6ps3pPmQL378e5bKdHtALRaOxMzAxERBtjg8qSCiZ8nrJl4cMPYfFiFaW8dKnFxbQbkU8iOTjvIIGTAnlwRSXgzuaZDd9+vlTtVhUX99QHRexFpVyPCST7FBVjkJrsMOtPreeD1R8QaYikZ7WezHhnBg4idUaLu3eVAoH0O/sArUA0GrsTu4QVZaCAeyEgZUtYMYwaBStWwLJl0L8/VIm/FlGGIexuGHumqXTqT0OeApC7TG78B/lTsX1FHF1S78/6BGXnmIrytiqJCgh8I5X9/Xr8Vzr81oFoYzS9a/ZmSqMpKfL4ep4JE+DRI1XOuHbtVHdjdbQC0WjsjLNz3AwkZgkrNQqkSBH4/HP15DpkCGzaZFExbUboJVM69blx6dQL1SyE/yB/yjQvk2g69eTwNyqu4xIqCLA/KhVJ6hbAYPHhxXRZ3wWjNDLYfzD/V///0qQ8rl6N86aLif9Ir2gFotHYGXMjumc25XJ69cHVVPX19dfKDrJ5M2zfDvUyUMjtzcM32TVhF8dWHkMalCtuqSal8B/kj1cdrzTdlAFCUPmqFpm2K6NmHdXS0OesfbP49I9PkUhGvTmKYa8PS7Ocw4dDRIRyjPDxSVNXVkcrEI3GzsQY0UWUkYJZC+Lk4MStJ7cIjw7HzSllKVdz51ZG9GHDVLW63bvT7/o5KFfcSzsuETA+gPObVbUGBycHKnasiN8AP/JVTDqdepJjoMqf9gLuAG7ASFQm3ZT7RsUxMWAiA7cOBGDCWxMY4D8gbYICR4/CokXg5KQ86tI7WoFoNHbG2cyI7ujgiGc2Ty6FXuLqg6uU8iiV4v6++krFD+zdq+pHtEppEIMNMBqMnFp3ioDxAVzfex0wueJ2N7nieqXOFfd5rqKy5m40bb8BzAFSflXjkFLyzT/fMPa/sQgEM96ZwSc+lknHN3iwKln7ySdQsqRFurQqWoFoNHYmJpWJiFKpN7yye3Ep9BKXH1xOlQLJkkVVrfv0U7Wk1bw5OKflUduCRIdHc3jJYXZN3MW9s/cAyJw7MzW+rEH1z6un2hX3eYzAz6jiTo+B7MD3KI+rtARzG6WRrzZ9xU97fsJROLKoxSI6VLJMdsMdO+DPP8HdXXnSZQS0AtFo7IyTmREdoEh2labvcujlVPf58ceqVsjZszB/vv2r14U/CGffz/vYPXU3j28+BiBH0Rz49velSpcqOGe2nIY7icrIusu0/R7wE1Awjf1GG6Pp/nt3Fh5aiIujC7+2+pXmZZqnsVeFlGrpEdTfvHkt0q3V0QpEo7Ez5ktYEKdArjy4kvo+ndUa+gcfqCDDDh3Uk62teXT9EUE/BLFv5j4iH0UCkM87H/6D/CnfujwOTpZL7hSJiuEYa3pfAJiGUiBp7tsQSYffOrD6xGoyO2dmfdv1vFU8wZp5KWbVKrXkmD+/SoqZUdAKRKOxMzFLWA4GiZSSIjlMM5AHqZ+BgLJ91KgBe/aouIJR8RZ9tg53T90lYGIAR5YcwRilFGOxesXwH+Sf4qy4ySEINes4btruDkxAZdBNK2FRYbz/6/tsOreJ7K7Z+bPDn/gV9rNAz4rISOV2DUrZZ8lisa6tjl0UiBAiF8oxoijKHfsDKeX9eNoZUClqAK5IKZuZ9hcDVqDS1OwHOkopI60vuUZjeZyEwOgocDBIjNHG2Cy8aVUgQqiYEH9/+P576N7d+okWg4OCCRgfwKn1p5T7k4ByrcrhN9CPQtULWXy8x8BQ1BKVRAUEzgHetFD/DyMe8u6yd/nvyn/kzpybvz/8myoFLBuhOWsWXLgAr72mlh4zEvZKDjwY2CalLAVsM23Hx1MpZWXTq5nZ/vHAFCllSVRt+wx22TWaOJwAo2kWYow2WsQGEoOfn1rGevpUGdStgZSSM3+cYeEbC5nnO49T607h6OJItZ7V+PLMl7Re1doqyuMvoDzwI+pGNhg4guWUx92wu9RfXJ//rvxHoayF+K/LfxZXHvfuxZWqHTdOue9mJOwlbnPivudFwA5gUHJONNVBrwe0Nzt/JMrpQqPJcDgCBmdHnCIMGKPiZiBXH17FYDSkOhlfDOPHw/r18Msv8OWXalnLEhgiDRxbeYxdE3Zx+9htAFyzu1L9s+rU7FUT9/zWMbrcRJWWXWnargrMQwUGWorgh8E0+qURJ+6coHjO4mzrtI2iOYpacATFyJFKidStq7zlMhr2UiD5pJQxuRpuAglFC7kJIfah0vSPk1KuQy1bhUopo01tggHLP95oNDbCkWdnIJmcM5E3S15uP7nNzcc3KZQtbT/vokVVbMj48cpA+99/aQsuDA8NZ//s/eyeuptH1x8BkLVQVmr1qUW17tVwzZbyuhfJwYhanhoEPEClHvkWpUwseSM7ffc0DX9pyJUHVyiXpxxbOm6hYNa0+nC9yMmTMGOGKhL1ww/pO+AzIaymQIQQW4m/gNdQ8w0ppRRCyAS6KSKlvCaEKA5sF0IcRf12UiJHD6AHgJdX0hXeNBpbY76EZTCLBbn95DaXH1xOswIBtXy1YAEEBKi0761bp7yP0EuhBE0N4uDcg0Q+VibHPOXz4NvPl0odKqUpuWFSHEPlr4pxzW2CKvRU1MLj7Lu+j8ZLG3M37C6+nr5sbL+RXJlyWXgURd++YDBAjx6qpktGxGoKREqZoI+bEOKWEKKAlPKGEKIAcDuBPq6Z/l4QQuwAqgBrgBxCCCfTLMQTuJaIHLNRpY3x8fFJSFFpNHbDETA6x81AQLny7ru+j8uhly3i8ZMtm0rM17OnijNo2hTckpkl5dreawROCuTEqhOx5WKL1S+GX38/SjQqYXGPKnOeAqOBiahliPwom0crkl+XPLlsvbCVlitb8jjyMY1LNmZV61VkcbGOS9Rff6lklzHfS0bFXkb0Daja9Jj+rn++gRAipxDC1fQ+N+APnJBSSuAf4qpMxnu+RpNReH4JC8yCCdPoiWVO165QoQJcugRTpybeVholp38/zcI3FjK3xlyOrzyOcBBU+rASPQ/2pNPWTpR8u6RVlccWoALwHWBA1eo4BbTG8spj1fFVNFnahMeRj+lQsQPr2663mvKIioqL9Rg+POMEDcaHvWwg44BfhRAfA5eBDwCEED7AJ1LKbkBZYJYQwohSdOOklCdM5w8CVgghxgAHUTY0jSZDEmNEB2JjJmJiQdISTPg8Tk4qOr1hQ1Wk6MMPodBzq2NRT6M4vPgwQVOCCDkdAoBrNleq9axGzV41yeaZmlJLKeM20AdYZtqugFpC8LXSeD/v/ZnP//wciaR3zd5MbjQ51YWgkjXez3DqFJQqpZwaMjJ2USBSyhDghSKNUsp9qHggpJS7gIoJnH8BsJAviUZjX5534wUsFgvyPA0aQMuWsHatKjq1fLna/+TOE/ZO38ve6XsJuxsGQHav7NT8qiZVP65qNcO4OUZgPjAQ5ZufCRhB2rPmJoSUktH/jmbEjhEAjK03liG1h1h1VhUSEue2O2kSuKS+mGK6IIN5HWs0Lx/mS1gxRnRLxoI8z5Qpav19xQro8PZd5K5Ajiw+QnS4cmws6FMQ336+lGtVzqKpRhLjJMpI/p9puyHKL99apd2N0kjvv3ozbe80HIQDM9+ZSfdq3a00WhxDhqg652+9Be++a/XhrI5WIBqNnYnXiG6WzkRKadGnYi8vydcdrnBi7i72f3Qmdn/ppqXx6+9nkeJNySUc+D/UmnYUkBf4AWiL5e0cMUQaIum8rjMrjq3AxdGF5e8v572ylsiYlTiBgTBnjspT9tNPGdNt93m0AtFo7Ex8S1g53XLi7uLO48jH3A+/bxFXUmO0kRNrThD4fSCGfdd5DYjGEbea3vRc6EvuMrnTPEZK2A58Apw1bXdHpZjIacUxH0Y85P1f32frha1kdcnK+rbrqVusrhVHVERHq/T6AAMGQJkyVh/SJmgFotHYmfiM6EIIvLJ7ceLOCS6FXkqTAol4FMHBeQcJ+iGIB5dVGFXm3JnJ2ag6Xy2tDsey0NmGmXpvAQOAJabtcsAsoLaVx7328BpNljXhyK0j5M2Slz/b/0m1gmkpaJt8fvoJDh9WQZ1DhybZPMOgFYhGY2fic+MFKJO7DCfunOD47eNULVA1xf0+DH7I7p92s3/WfiIeRACQq1QufPv54t3JG+dMzmwKV1UL+/aFX3+1yMdJEANKUXyNigZ2BYahjObWtiUfv32cxksbc/XhVUp7lGZTh00Uy1nMyqMqgoOVuy4oRZLZMjWz0gVagWg0dsaBOAUSbaZAvPN589vJ3zh86zAd6Zjs/m4evkngpECOLT8W59VVxwu//n6Ufrc0wiFu8X3yZBXUtmoVbNmivLSswT5UHMc+0/bbqFodJawz3DP879L/aL6iOQ8iHuBX2I8NbTfgkdnDBiMr+vSBx4+hRYuXw3BujlYgGo2dEcQZ0aNNXligFAjA4VuHk+xDSsn5zecJnBTIha0XVL8OgvIflMe3ny+FasSfDsXLS5VPHTJE1eE+csSy9ShCUbmLfkalW/cEpgItsZ6R3JyVx1bSaV0nIg2RtCzTkqXvLSWTcyYbjKzYtEmljsmcOengzYyIViAaTTpAxjMDqZy/MgCHbh5K0BMrOiKao8uOEjQ5KDYjrnMWZ6p2q0rN3jXJWSxpk3S/fioe5MgRpUwmT7bA5wF+AfqjAgMdUcGBIwBbmFuklEwOnEz/Lf0B+LLGl0xpNCXNmY1TwtOn8MUX6v3IkUpZv2xoBaLRpAdMCiQyKk6BeGX3IodbDu6G3eXG4xvPZIR9eu8p+2buY89Pe2JrjGctmJUavWpQrUc1MuVM/lO2s7Oqm16jhsoK+8EHUKtW6j/KCeAz4H+m7drADBKICrYCBqOBvpv78uOeHwGY2GAi/Xz72cw1OYZhw+D8eZU+5quvbDq0zdAKRKNJBxhNXlgGsxmIEIJK+Srx7+V/OXzzMAWzFuTe+XsE/RDEofmHiAqLAiBfpXz49vOlQtsKqc6IW62aikyfMEFVxTtwAFxTGHz+BJX4cBIq8WFuVBLEzthmuQrgadRTOq7tyJqTa3B2cGZRi0W0q9jORqPHsWuXCth0cFDK2dkaofTpAK1ANJr0QDxLWKDsIP9e/peD2w7yaMgjTq09FZsRt0SjEvj286X4W5apMT5ypEpxcuIE/N//wbffJu88icpm2hu4glIWPVEBgtZJhB4/d8Pu0mJFCwKuBpDdNTtr26y1SYzH8zx9Cl26gJQwaBBUr25zEWyGViAaTTogPiO60WCk9PHSdJ3XlairUZzkJA7ODnh38qZW31rkq5hQHbbUkSkTzJ0Lb7yhFMj77yddp+Ii0AvYaNqugjKY17SoZElz8s5J3l3+LhfuX8Azmyd/dfiLCnkr2FgKxfDhcOYMlCsXl/fqZUUrEI0mPWA2A4l8EsmhhYcImhLE/fP38cKLyEyR1PuqHjW+qEHWglmtJsbrr8Nnn6lKeV27QlBQ/HW6n6Kixsej0pFkA8agXHVtfVPZemErrX5txYOIB1QrUI0N7TZYpYJgcggMVEkSHRxUAa+ULgNmNLQC0WjSAyYFcmrBIfb038LTe08ByF4sOyvLrORglYMMGznMJi6o48bB77/D/v3q/bBhcccksA7lURWT5rE98D1QwOqSvcjs/bP57I/PMEgDLcu0ZEnLJVar45EU5ktXAwdarvZ8esZeBaU0Go0Z0rSEdeu/Kzy995RCNQvRelVrep3txcN3HxLhEsGx28dsIkvWrDDPVGFn5Eg1CwFVzKkR8B5KeVRCeVotxfbKw2A00G9zP3pu7IlBGhjkP4jVH6y2m/IAtXR1+jSULfvyL13FoBWIRpMOuFejEEYHQeEWZeiyswsfB36s0qk7OqQooNBSNGig4kMMBmjXE76MUG64W4AcqCjy/cDrNpMojseRj3nv1/eYHDQZJwcn5jWbx7i3xlm1CFRSbNv27NJVcssFZ3TssoQlhMgFrASKApeAD6SU959rUxeYYrarDNBWSrlOCLEQeAOVUgfgIynlIetKrdFYj4tdqrC1kzenHB14Pt7MO583y48t5+CNgzaVacxYWOUGlz6Haa7Ku6o7MBbIY1NJ4rjy4ArNVzTn0M1D5HTLyW9tfuPNom/aSRrFrVuquqOUMGIE1LS1B4EdsZfKHgxsk1KWAraZtp9BSvmPlLKylLIyUA8IA/42azIg5rhWHpqMjiMgHR0wxHPMr7AfAP9c+sdm8hwA6rvClTGo9alAGL1ZlZa1l/L49/K/+Mz24dDNQ5TMVZKgbkF2Vx5GI3TuDDdvKu81c3vRq4C9FEhzYJHp/SKgRRLtWwF/SSnDrCmURmMvYpYCouM5VsuzFu4u7py8e5Lgh8FWlSMEVaPDB9gF5AM6/wP4w4QP4NIlqw4fL1JKZuydQf3F9bkTdocGxRuwu9tuSnuUtr0wzzFpEmzeDB4esHQpONouU0q6wF4KJJ+U8obp/U3U7zQx2gLLn9s3VghxRAgxRQiRoLOcEKKHEGKfEGLfnTt30iCyRmM9Yu478c1AnB2dqVtUBcRtOb/FKuNHAT8CpVAp1x1QnlangQVvQvNm8PChWqqJjk/LWYmI6Ah6/N6Dz//8nGhjNP18+/Fnhz8tUmArrezeDV9/rd4vWgSF4s9X+VJjNQUihNgqhDgWz6u5eTsppUR5BybUTwGU/W6z2e4hKJtIdVSw66CEzpdSzpZS+kgpffLksdfkW6NJnMQUCECD4irP+t8X/k6gReqQwB+of7DewH2gPnAEmAxkR5VenTsXChaEgIC4m6a1ufHoBnUX1WXuwbm4ObnxS8tf+L7h9zg52D/6IDQU2rZVyrRPH3jnHXtLZB+s9k1IKd9K6JgQ4pYQooCU8oZJQdxOpKsPgLVSyiizvmNmLxFCiAWopJ8aTYYlsSUsgIYlGgIqaM4ojRbxODoO9CXOsFgKFc/RlBdzV+XODStWQN26MHGiMhS//36aRUiQPdf20HJlS64/uk7hbIVZ22atzaoHJoWUKl/YpUsqh9h339lbIvthryWsDagca5j+rk+kbTueW74yKR2ESgDUArCNg7xGYyWSmoGU9ihN4WyFuRt2l0M3D6VprDuobLmVUMojO2q2cQxoRsKJD+vUUcoD4KOP4NSpNIkRL1JK5uyfw+sLXuf6o+vU8arDvh770o3yAJXm5bffIFs2pVRf9mjzxLCXAhkHNBBCnAXeMm0jhPARQsyNaSSEKAoUJi4zdAxLhRBHgaOopJ9jbCG0RmMtklIgQojYWUhq7SCRKEVRCpWvSgCfA+dQ9o7klJX96iuV7v3xY2jeHO7fT/KUZPMk8gmd13Wmx8YeRBgi+NTnU7Z22kreLHktN0ga2bhR1UwRApYtg5Il7S2RfbGLApFShkgp60spS0kp35JS3jPt3yel7GbW7pKUspCU0vjc+fWklBWllBWklB9KKR/b+jNoNJYkZgkrIQUCqbeDSOA3oDzQDxU81Qhl55iGegJLLkKoKHVvb5Uw8IMPLGNUP3HnBDXm1mDJkSVkds7M4haLmfHODFwcrV0tPfkcPw4dOqglrDFjXl27hzk6El2jSQfEzEASuxfXL14fgWDnlZ3cf5q8R///AD/gfdRMoyzwJ7AJKJdKWd3dYcMGyJsXtm5VsxKZoBtM0vxy5Beqz6nOiTsnKJu7LHu67aGjd/JrwNuCmzeVwnj4EFq1UiWANVqBaDTpgqSWsAByZ85NvWL1iDREsuLYikT7O4GyZ7wOBAF5genAYaBx2sXFywvWrQMXF5g+Hb7/PuV9PI16So/fe9BxbUfCosLoULEDe7rvoXze8haQ0HKEhUGzZnD5snIeWLxYzcQ0WoFoNOmC5CgQgI+rfAzAvIPz4j1+DZVupCLwO5AFVYf8HMpwbsnCeL6+6mYKKvvs0qXJP/dsyFn85vsx58AcXB1dmf3ubJa0XIK7iy0qpiefqCho0wb27oWiRWH9elU3RaPQCkSjSQck5cYbQ4syLcjhloP9N/Zz+GZccsUHwNcoA/lclIH8U5TiGAlYq4JImzYwebJ636UL/PVX4u2llMw9MJcqs6pw6OYhSuQsQVC3ILpX627zmuVJYTSqz7RxI+TKBX/+CfksW8Mrw6MViEaTDkjuDCSTcybaV2gPwIJDCwgHfgBKAN+hCj29j1rCmgHkt4awz9Gnj8rcGxUF770H27fH3+5u2F3e+/U9uv/enSdRT2hboS37e+yncv7KNpAyZUgJvXqpWZW7u1KMZcvaW6r0h1YgGk06ILkKBKBrla7g4MxcB2dKSkkfVA6rOkAgsBqwdZaoiRPhk08gPByaNoX//nv2+OZzm6n4c0XWnVpHNtds/NLyF5a/v5zsbtltLGnSxCiP6dOVjWf9+lejOFRqsH9OAI1Gk+wlrCjgYIGqOH91iSdZC/IEFRA4BniXhIMArY0Q6oYbHg4LF8Lbbysje83XHzJwy0Bm7Z8FQB2vOixuuZiiOYraSdLEMRrhyy9VSV9XV1i7FurVs7dU6RetQDSadEBSMxADsAz4FjgvBGQtCHdOkHfvdHY3moJbOoiXcHBQObNAKZEmvf4ie8cehEQF4+zgzMg3RzLIfxCODukzZW1UlKoD/8svSnmsW6cUoSZhtALRaNIBCSmQaGAFqohTTOaQUsA3RgNjV7Xm9J0TzPEow5c1v7SRpInj6AgTp90jMH8fTrstJiQKijhVZ2P3+VTIW8He4iXI48cqvmPzZsiSRc08GjSwt1TpH20D0WjSAc8vYYUDM1G2jI4o5VEMWIAykHd0cGR8vf8DYNS/o3gY8dCm8saHURpZfHgxFX4ux2m3xThJN/h7ApeH72LB+AoYkmPgsQOXL6s8X5s3Q5488M8/WnkkF61ANJp0QMwM5AEwCSiOcsO9iFIi81FK5CPilE2z15rhX9ifu2F3GbdznG0Ffo4DNw5Qe35tOq/rzK0nt6jjVYfjXx5mfrcBODk4MXmyiuS+e9euYr7Ajh3g4wOHDqm8VgEBUL26vaXKOGgFotGkA2IUSF9UbYIbgDewEjXj6MKLyQ6FEExsoNLjjg8Yz/aLCfjPWpG7YXfp+XtPfGb7EBgcSL4s+VjUYhE7PtpBaY/SdOkCf/+tKvZt3gxVqkBgoM3FfIHoaBg1CurXV0qtUSPYswdKlbK3ZBkLrUA0mnSAm9l7P1SRp4OoYjiJmZx9C/sytM5QjNJI29VtufrgqjXFjOVB+ANG/W8UJX8syewDs3F0cKSfbz/OfHmGTt6dnqlXUrcuHDyoIteDg6F2bZVLKiLCJqK+wPnzSqYRI5TX1eDB8McfkDOnfeTJ0EgpX5lXtWrVpEaTHjkppfxcSrlDSmlM4bnRhmjZYHEDyUhkjTk15KOIR5YX0MTD8Idy7L9jZc5xOSUjkYxENlzSUJ64fSLJcyMipBwwQEohpAQpy5WT8p9/rCZqvOOPGSOlm5sav0ABKbdssd34GRlgn4znnipkWtJoZjB8fHzkvn377C2GRmNx7obdpdrsalx5cIUq+auwsf1GCmYtaLH+Q8JCmLV/FlOCpnA3TBkyXi/yOqPeHMUbRd9IUV+7dqkUIWfOqO2WLWHcOChtpehHoxF+/RWGD4ezZ9W+Dz+EKVNUpUVN0ggh9kspfV7YrxWIRvNycDbkLE2WNeHcvXN4ZvNkxfsr8PfyT3V/UkqCgoOYd3AeS48uJTw6HABfT19G1x1NvWL1Up2/6ulTlUPru+/gyRMViNiqFQwapMrEWoLwcFi5Uo1z5IjaV7q0ChKsX98yY7wqpCsFIoRojcrxVhaoIaWM964uhHgbmIpaBp4rpYypXFgM5R7vAewHOkopI5MaVysQzcvO3bC7NF/RnF1XdwHQqlwrRtcdTZncZZJ1vlEaOXDjAOtPrWfViVWcDjkde6xxycZ8VesrGhRvYLHEh9evK1vEokUqkA+gcmVVMrd5c5UBNyUYDGqGs3q1qhgY4/Xl6anG6dwZnC2ZkvgVIb0pkLKAEZgF9I9PgQghHIEzQAMgGNgLtJNSnhBC/Ar8JqVcIYSYCRyWUv6c1LhagWheBcKjwxn771gmBU7iafRTALzzedPstWZUzFuRErlK4O7ijlEaeRTxiEuhlzh77yy7r+0m8Gogd8LuxPaV3z0/HSt15OMqH/Na7tesJnNwsJopLFz4bJnc0qXh9dehYkUoV07FaXh4qBlLZCSEhMDFi2o5LChIeXiFhMSdX6WKSk3Srh24ub0wrCaZpCsFEju4EDtIWIH4AiOllI1M2zE1wMYBd4D8Usro59slhlYgmleJ4IfBfLvjW1YeX8mjyEfJPs8zmyfNSjej2WvNqF+8Pk4OtktYERGhqh0uXw7btqkKgCmlWDG1HNaqlYrpSGdZ4jMkCSmQ9JzKpBBg7pMYDNRELVuFSimjzfYXSqgTIUQPoAeAl5eXdSTVaNIhntk8mdNsDtOaTGPrha38c+kfzt07x/n754mIjsBBOJDJORNFcxSlWI5iVC1QFV9PX4rnLG632hyurtC6tXpFRanYjP374ehRZQAPCYmbYbi6QvbsapmrWDGlLPz8oEgRrTRshdUUiBBiK/GXIxgqpVxvrXGfR0o5G5gNagZiq3E1mvSCq5Mr75R+h3dKv2NvUVKEszP4+6uXJn1iNQUipXwrjV1cAwqbbXua9oUAOYQQTqZZSMx+jUaj0diQ9ByJvhcoJYQoJoRwAdoCG0xBLf8ArUztOgM2m9FoNBqNRmEXBSKEaCmECAZ8gT+EEJtN+wsKIf4EMM0uvgA2AyeBX6WUx01dDAL6CiHOoWwi82z9GTQajeZVRwcSajQajSZREvLCSs9LWBqNRqNJx2gFotFoNJpUoRWIRqPRaFKFViAajUajSRWvlBFdCHEHuJzK03MD6awgJ6DlSilarpSh5UoZL6tcRaSUeZ7f+UopkLQghNgXnxeCvdFypQwtV8rQcqWMV00uvYSl0Wg0mlShFYhGo9FoUoVWIMlntr0FSAAtV8rQcqUMLVfKeKXk0jYQjUaj0aQKPQPRaDQaTarQCkSj0Wg0qUIrEDOEEK2FEMeFEEYhRIIub0KIt4UQp4UQ54QQg832FxNC7DbtX2lKQ28JuXIJIbYIIc6a/uaMp01dIcQhs1e4EKKF6dhCIcRFs2OVbSWXqZ3BbOwNZvvteb0qCyECTd/3ESFEG7NjFr1eCf1ezI67mj7/OdP1KGp2bIhp/2khRJJlmy0sV18hxAnT9dkmhChidize79RGcn0khLhjNn43s2OdTd/7WSFEZxvLNcVMpjNCiFCzY1a5XkKI+UKI20KIYwkcF0KIH00yHxFCVDU7lvZrJaXUL9MLKAu8BuwAfBJo4wicB4oDLsBhoJzp2K9AW9P7mcCnFpJrAjDY9H4wMD6J9rmAe0Bm0/ZCoJUVrley5AIeJ7DfbtcLKA2UMr0vCNwAclj6eiX2ezFr8xkw0/S+LbDS9L6cqb0rUMzUj6MN5apr9hv6NEauxL5TG8n1ETAtnnNzARdMf3Oa3ue0lVzPtf8SmG+D6/U6UBU4lsDxJsBfgABqAbstea30DMQMKeVJKeXpJJrVAM5JKS9IKSOBFUBzIYQA6gGrTe0WAS0sJFpzU3/J7bcV8JeUMsxC4ydESuWKxd7XS0p5Rkp51vT+OnAbeCHS1gLE+3tJRN7VQH3T9WkOrJBSRkgpLwLnTP3ZRC4p5T9mv6EgVPVPa5Oc65UQjYAtUsp7Usr7wBbgbTvJ1Q5YbqGxE0RK+S/qYTEhmgOLpSIIVc21ABa6VlqBpJxCwFWz7WDTPg8gVKpCWOb7LUE+KeUN0/ubQL4k2rflxR/vWNMUdooQwtXGcrkJIfYJIYJiltVIR9dLCFED9VR53my3pa5XQr+XeNuYrscD1PVJzrnWlMucj1FPsjHE953aUq73Td/PaiFETOnrdHG9TEt9xYDtZrutdb2SIiG5LXKtrFYTPb0ihNgK5I/n0FAppd1K4yYml/mGlFIKIRL0vTY9XVREVXKMYQjqRuqC8gcfBIyyoVxFpJTXhBDFge1CiKOom2SqsfD1WgJ0llIaTbtTfb1eRoQQHwI+wBtmu1/4TqWU5+PvweL8DiyXUkYIIXqiZm/1bDR2cmgLrJZSGsz22fN6WY1XToFIKd9KYxfXgMJm256mfSGo6aGT6SkyZn+a5RJC3BJCFJBS3jDd8G4n0tUHwFopZZRZ3zFP4xFCiAVAf1vKJaW8Zvp7QQixA6gCrMHO10sIkQ34A/XwEGTWd6qvVzwk9HuJr02wEMIJyI76PSXnXGvKhRDiLZRSfkNKGRGzP4Hv1BI3xCTlklKGmG3ORdm8Ys5987lzd1hApmTJZUZb4HPzHVa8XkmRkNwWuVZ6CSvl7AVKCeVB5IL6sWyQyjL1D8r+ANAZsNSMZoOpv+T0+8Laq+kmGmN3aAHE67FhDbmEEDljloCEELkBf+CEva+X6btbi1ofXv3cMUter3h/L4nI2wrYbro+G4C2QnlpFQNKAXvSIEuK5BJCVAFmAc2klLfN9sf7ndpQrgJmm82Ak6b3m4GGJvlyAg15diZuVblMspVBGaUDzfZZ83olxQagk8kbqxbwwPSAZJlrZQ3PgIz6Alqi1gIjgFvAZtP+gsCfZu2aAGdQTxBDzfYXR/2DnwNWAa4WkssD2AacBbYCuUz7fYC5Zu2Kop4sHJ47fztwFHUj/AVwt5VcgJ9p7MOmvx+nh+sFfAhEAYfMXpWtcb3i+72glsSamd67mT7/OdP1KG527lDTeaeBxhb+vScl11bT/0HM9dmQ1HdqI7m+A46bxv8HKGN2blfTdTwHdLGlXKbtkcC4586z2vVCPSzeMP2Wg1G2qk+AT0zHBTDdJPNRzLxLLXGtdCoTjUaj0aQKvYSl0Wg0mlShFYhGo9FoUoVWIBqNRqNJFVqBaDQajSZVaAWi0Wg0mlShFYhGo9FoUoVWIBqNRqNJFVqBaDR2RAhR3ZQU0E0IkUWo+iQV7C2XRpMcdCChRmNnhBBjUNHomYBgKeV3dhZJo0kWWoFoNHbGlFtpLxAO+Mlns7hqNOkWvYSl0dgfD8AdyIqaiWg0GQI9A9Fo7IxQNbJXoIoQFZBSfmFnkTSaZPHK1QPRaNITQohOQJSUcpkQwhHYJYSoJ6XcntS5Go290TMQjUaj0aQKbQPRaDQaTarQCkSj0Wg0qUIrEI1Go9GkCq1ANBqNRpMqtALRaDQaTarQCkSj0Wg0qUIrEI1Go9Gkiv8HuJT+9B1gRK0AAAAASUVORK5CYII=" }, "metadata": { "needs_background": "light" } } ], "metadata": {} }, { "cell_type": "markdown", "source": [ "This nicely shows the shock developing in the center of our domain, which forms from the collision of the two initial velocity \"bumps\", the positive one on left (moving right) and the negative one right of the center (moving left).\n", "\n", "As these lines can overlap quite a bit we'll also use a different visualization in the following chapters that shows the evolution over the course of all time steps in a 2D image. Our 1D domain will be shown along the Y-axis, and each point along X will represent one time step.\n", "\n", "The code below converts our collection of velocity states into a 2D array, repeating individual time steps 8 times to make the image a bit wider. This is purely optional, of course, but makes it easier to see what's happening in our Burgers simulation." ], "metadata": {} }, { "cell_type": "code", "execution_count": 6, "source": [ "def show_state(a, title):\n", " # we only have 33 time steps, blow up by a factor of 2^4 to make it easier to see\n", " # (could also be done with more evaluations of network)\n", " a=np.expand_dims(a, axis=2)\n", " for i in range(4):\n", " a = np.concatenate( [a,a] , axis=2)\n", "\n", " a = np.reshape( a, [a.shape[0],a.shape[1]*a.shape[2]] )\n", " #print(\"Resulting image size\" +format(a.shape))\n", "\n", " fig, axes = pylab.subplots(1, 1, figsize=(16, 5))\n", " im = axes.imshow(a, origin='upper', cmap='inferno')\n", " pylab.colorbar(im) ; pylab.xlabel('time'); pylab.ylabel('x'); pylab.title(title)\n", " \n", "vels_img = np.asarray( np.concatenate(vels, axis=-1), dtype=np.float32 ) \n", "\n", "# save for comparison with reconstructions later on\n", "import os; os.makedirs(\"./temp\",exist_ok=True)\n", "np.savez_compressed(\"./temp/burgers-groundtruth-solution.npz\", np.reshape(vels_img,[N,STEPS+1])) # remove batch & channel dimension\n", "\n", "show_state(vels_img, \"Velocity\")" ], "outputs": [ { "output_type": "display_data", "data": { "text/plain": [ "
" ], "image/png": "iVBORw0KGgoAAAANSUhEUgAAA2AAAAEeCAYAAADsJXY8AAAAOXRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjMuMywgaHR0cHM6Ly9tYXRwbG90bGliLm9yZy/Il7ecAAAACXBIWXMAAAsTAAALEwEAmpwYAABfaElEQVR4nO29e7hud1Xf+xm/Od+11w63ANEQkiho06NUEWyKtpwq5dLGyyE+pxiDgmBDc2xF7U0J4gOUyjnR9tHSowfdD6DBIiHihdhSEbkcnx4FkyiWm0oIlyQNxEDue6/LnL9x/phzJcvNvmTv/fuutcZ+xyfPerLed71rrN+ea653/r5zjPEd5u4kSZIkSZIkSZIkespuLyBJkiRJkiRJkmRZSAGWJEmSJEmSJEmyQ6QAS5IkSZIkSZIk2SFSgCVJkiRJkiRJkuwQKcCSJEmSJEmSJEl2iBRgSZIkSZIkSZIkO0S/2wtIkiRJkiRJkiT5Rxc92b9wx30n/H033PCpd7n7RYIlSUgBliRJkiRJkiTJrvOFO+7jg9f/uxP+vt5ecJZgOTJSgCVJkiRJkiRJsus4Tq3jbi9DTgqwJEmSJEmSJEn2AI77sNuLkJMCLEmSJEmSJEmS3cfBPTNgSZIkSZIkSZIkchynZgYsSZIkSZIkSZJkJ8gSxCRJkiRJkiRJkh0iBViSJEmSJEmSJMnO4I7XFGBJkiRJkiRJkiQ7Q2bAkiRJkiRJkiRJdoIsQUySJEmSJEmSJNkhHOrmbi9CTgqwJEmSJEmSJEl2HffMgCVJkiRJkiRJkuwQDmnCkSRJkiRJkiRJsgN4CrAkSZIkSZIkSZKdI0sQkyRJkiRJkiRJ9BiOZQYsSZIkSZIkSZJkB1iSEsSy2wtIkiRJkiRJkiRZFjIDliRJkiRJkiTJHmA5MmApwJIkSZIkSZIk2QM4liYcSZIkSZIkSZIkO4ADddztVchJAZYkSZIkSZIkyR4gXRCTJEmSJEmSJEl2CM8MWJIkSZIkSZIkyY6wJDb0KcCSJEmSJEmSJNkTWGbAkiRJkiRJkiRJdgDPEsQkSZIkSZIkSZIdIzNgSZIkSZIkSZIkO0JmwJIkSZIkSZIkSXYEc1+KDFjZ7QUkSZIkSZIkSZIAUwbsRD+Og5ldZGZ/YWY3mtkVR/j6z5nZh+aPvzSzu7Z9bdz2tWtb/BMzA5YkSZIkSZIkye4jyICZWQf8AvAc4BbgOjO71t0/9uCP9X+57fU/DDx1W4hD7v6UlmvKDFiSJEmSJEmSJHuD9hmwpwE3uvtN7r4BXA1cfIzXPx94a6N/zRFJAZYkSZIkSZIkyR7AsVpP+OM4nAvcvO3xLfNzX4KZfSXwROC9255eNbPrzewDZvZdp/CPe4AsQUySJEmOi5k9AfgUsHD34RTi3Ac82d1varW2JEmS5DTBOVkXxLPM7Pptjw+4+4GTiHMp8HZ3376Ir3T3W83sq4D3mtmH3f2TJ7PILVKAJUmSLBFm9rvAH7v7Kw97/mLgl4DzTkVgHQ93f/i2n/krwC3u/pOqn5ckSZIsBXe4+4VH+dqtwPnbHp83P3ckLgV+aPsT7n7r/P+bzOz9TP1hpyTAsgQxSZJkubgKeIGZ2WHPvxB4i1J8JUmSJMmxcUUP2HXABWb2RDNbYRJZX+JmaGZfAzwa+KNtzz3azPbNn58FPB342OHfe6KkAEuSJFkufht4LPD3t54ws0cD3wm82cyuMLNPmtkXzOwaM3vMkYKY2ePN7Foz++Js6/tPt32tM7OfmOPca2Y3mNn589fczP6GmV0OfB/w42Z2n5n9jpn9mJn9xmE/5z+Z2euaH4UkSZJkT2JeT/jjWMw3Fl8KvAv4OHCNu3/UzF5jZs/d9tJLgavd3bc997XA9Wb2Z8D7gCu3uyeeLFmCmCRJskS4+yEzuwb4fuAP5qcvAf4ceAbwXcC3An8F/Ccm697nHyHU1cBHgMcDXwO828w+6e7vBf7V/D3fDvwl8GTg4GHrOGBmf49tJYhmdg7wajM7093vMrOe6YL4bW3+9UmSJMmexv1ke8COE9bfCbzzsOdeedjjVx/h+/4Q+PrW68kMWJIkyfJxFfA8M1udH3///NwPAq9w91vcfR149fy6v3azbs5mPR14mbuvufuHgDfMcQBeAvyku/+FT/yZu3/heIty99uYROF3z09dxFTXf8Mp/FuTJEmSSNR64h/BSAGWJEmyZLj7fwfuAL7LzL6aaUbKrwFfCfyWmd1lZncxlWqMwNmHhXg88EV3v3fbc5/hQVvf8zn5BuWrgBfMn78A+NWTjJMkSZJEwz0FWJIkSXLa8mamjNULgHe5++eZ5qR8m7ufue1jdcsBahv/E3iMmT1i23NfwYOuUjcDX/0Q1uBHeO63gSeb2dcx9aW95SH/i5IkSZLwWB1P+CMaKcCSJEmWkzcDzwb+KVPWCeAXgdfOgygxsy+b7en/Gu5+M/CHwP9lZqtm9mTgMuA/zy95A/DvzOwCm3iymT32CGv4PPBVh8VeA97OlJH7Y3f/7Kn+Q5MkSZIoZAYsSZIkOU1x908ziaiH8aAd7+vmz3/PzO4FPgB801FCPB94AlM27LeAV7n7789f+1ngGuD3gHuANwL7jxDjjcCT5pLH3972/FVMTc9ZfpgkSbJMOEshwOyvOy0mSZIkye5iZl/B5Mr4OHe/Z7fXkyRJkuwMf/tvrvgH/+/HnfD3LS66+YZjDGLec6QNfZIkSbJnMLPCZGN/dYqvJEmSJcOR2NDvNfacADOzi5jKYDrgDe5+5S4vKUmSJNkBzOxhTH1hn2GyoE+SJEmWCMOxgCWFJ8qeEmBm1jEN/XwOcAtwnZld22LidJIkSbK3cff7gYfv9jqSJEmSXSQF2I7zNOBGd78JwMyuBi4GUoAlSZIkSZIkyenMlgnHac5eE2DnMs2P2eIWju7AxZmLff64fWc0X0TakiRJkiRJkiSR+dz6Qe7eXLfdXseJ4SnA9iJmdjlwOcDZ+/bzpm/8+7u8ooeGe5zzvwZaK8Q6ttGIdi4ksci/3SRJkgdp/Z74zz7ynqbxdgQH6umfCtlrAuxW4Pxtj8+bn3sAdz8AHAB40iMf5fv6zeaLUGwKVBtZ1QbGibMxiiYSIm06O1HcaL+zKEQ6twCw0/8iezqRf7dJosXyPXEiM2A7znXABWb2RCbhdSnwvUd7sZmzum+96QKqm2QTo9oYjaNmiywTdkEu4MqNRpRjADohrhJ2CnLTGeucjUak86vLzWEo8u82yTNg77KnBJi7D2b2UuBdTHu0N7n7R4/2ejNnsRiarqFWnQBTxC2l/V0C1Vqn2KV5zFoFx5VYFy/VWiNtDlXHINKmU5a5DnQMIp2zEOz8CnZsl53TP4eQHJ847y8Pkj1gu4K7vxN450N5rZnTL9qWIHo1am0vEtyLRCiYYK2gLMNs/0dVSpyMpUKAgkaEAligDVekzaFKJJjoYhvp2KoETaSybBmCYxtNMEdCc7VJQhHxzyt7wPY+Zs7KykbTmO4mKeurteCCTbJCLIJGKIxjHLGoyyhpMpZFcGhThOqIJGwhlgCLJm4VRPp9hcoARtvNBjq2KcSTv4Zgr7TXiC3ASqXf116AdQoBNhbJRXEYNL/CKhBLqkxVpA19CdZjqLkoat5YFecXxBOiClLc6pCcX4E23hBn8x1JhEMKcRWRhHiclW7HMwO21zFz+pXGJYhujJvtf/FmJdSGa7T2IlTVW1ZHzR+qmUIoac4BVXZRtZmNsumsbhJXKnfDbJTEVdAJXFN074ea9wNZma/o/IqEYvMd7RgoOP1zCA+N1udCNCEejixB3PuYOd2+xgJs6CQXxHHo8SBNhS7adI6icklTiUXBhss91puKatMZJQvYmYcRixBr3IVCgEKK0Dly84iRsqDJRJxbvjpUN9HiEPTfngJsb2PF6VfbliDWoWCCsr7SVUlZ31gEdw5FTpBl7DQZMNGbq6IXULU5NFEWcNkNTiIJUAiWsRSRIlR0M0I4O0IjRONs4CJlVyHW+0Go0kbFe0zziHrcl6IFLLYAw5zS2ITDFoWy2f43Xzc7TLChN5EAq4o+OFE2QZUBU5UGKTaIqrVKMlUyUdP+GJSiytbF6VuEWHMBZS6Igfoso228I5XjKoiVXYXMsAYi6mHNDNgepzildQni2EnqpguEaYz2ahJhZ6XggjJExT1kYzoXWuNumKK0UdT/pPh9uaC3DkRCSXZnepRlKaJkE2rV3eSIRKj+p7ShD4PqPUZFUVj5BiNIl4oeZykaGEMLMDOFANOUhFSg64K8GbpJRE1BIxQUuJtkrT4WSWF+rYLsatUcA1mJWG3vNGqdJmvrQvERJZuQYjGeWFQg23YH67NcdqIJRgUpQreRAmyPU5yyv7UAE2V/ulGSUVHg1SQXLx+rpgRR4FoJGiv+WkTngKIHrGos/ovK4CTQhkshFkFYKhiobFSBwhEV4pWiKohU3gqZsQThrL1Ix1aA4uZZWAOSoMs+EWILMANbbfxbGhxoK+pgzv7Uxn9ctYBiI6vabBWXxW6NV5PcmTWVq54AN5u6YRtj1Zrf9nbXCBpAVjKqEnYKFGWj7iaxc1a8x6j+bpWlqK3RiUVJ2FCb+UjXBQjmrhioiiOZcdX7uF0EvA7ogDe4+5WHff3FwL8Hbp2f+nl3f8P8tRcBPzk//1PuftWprie8AGPRuOSouOaKUAcYGl8VuyrZHKqc1CpTH1hzXFByZQVXGJwMGidIBW6iLI21N/cwRo1Qck3JKFV0J1212RAIO9mmU7Q7VInFKEIhXt+iAg/jWBjlOqNG4ugcNasUicbbcJvc2n4BeA5wC3CdmV3r7h877KVvc/eXHva9jwFeBVzIlJu7Yf7eO09lTbEFWDFYaSxqhhFTuK/UCr3ApU3UnyKbTqMoaWseEbzWqV+rMZX2pkSy35dodIJS2LUP6qIMmCibIPh9deYyYafYIsqygIHEYrT+iSh9ixCpTyeOWARtdrE1kcRt2KHR7c/dpwE3uvtNAGZ2NXAxcLgAOxL/CHi3u39x/t53AxcBbz2VBcUXYGesto05DFDaWtsDWKlzeWNjavtySasmyf5gDoI3LokLYjVc0a+luHPmRhWNI9D0O1RJeYFidMIUV5G11VwUO2u/81aVdyqELcTrBVQQyWkUdGWjrYlkcgOxSlGVgrE10UpGw6EpQTwXuHnb41uAbzrC6/6xmX0L8JfAv3T3m4/yveee6oJiCzAzWKy0jVmKZv5A3ZCYezAINlxjkZRhlpUBb12GCZTF0DwmVSNqihWIUm7TiURNNUyUAWqNuya7KJspJVitrLxTdA64yORFk2UeNf0kImEnEzVdjCyFoeljjrShV7oVRhKMYSzjY5xWh2EnmwE7y8yu3/b4gLsfOIHv/x3gre6+bmb/B3AV8MyTWchDIbgAK7DSWIANGvHBMEJpHLc6tuLNJ4Zb63XO+IbmAi65k140G+8KEoEv2SBXo3Tt7/ZWOijt46qEkmQcgWjGWukEN2SqSe4ZmGwHU5pv6JWGIaF6ARWIxKKMLBsNw5Sx1BzYOKWoS8Ud7n7hUb52K3D+tsfn8aDZBgDu/oVtD98A/My2733GYd/7/lNZKEQXYKXgq2e0jTlsajYGq3USYa2pG+1vcIzCKXitw1aT9cEpmGahxbg779UkoxOKqBRVcsaKsnVKYdc8ZqSePZDUJE+GoIK/skCb+U6VpcksYKgsIBBGLCozi1FKUcP2gLU/FtcBF5jZE5kE1aXA925/gZmd4+63zQ+fC3x8/vxdwP9pZo+eH/9D4OWnuqDYAswMX9nXNmbRZMAMoFeUytX2GZVaMck7bA3zxi27K2sucZms7dsWsQ6NE+TYBSoRi1WGqRFKomHUKiMSSVR0IxkEKMs7m8cMlAUETSbQ0PQDpljUlY0mQgQ9YO4+mNlLmcRUB7zJ3T9qZq8Brnf3a4EfMbPnAgPwReDF8/d+0cz+HZOIA3jNliHHqRBcgBV8pbEJRymaDNjQ3iwDgF6wORxGTZkcjgs0qGIIsbWe2bYdxYW2V4xOMDQWJ0wphQCoyjBVYw6aj7oArFMJME1cWc+eQtmJjq2yvDNCSECmxGXLFZX+NydYyai7afrvFShOgRiX2i9FcJ65+zuBdx723Cu3ff5yjpLZcvc3AW9quZ4dF2Bmdj7wZuBsJj/9A+7+utln/23AE4BPA5ccz2PfzaiNBZiVDleIpZV9mhlY49A+YyczIhnbvxFW5ixg27CTQFCUjBoILgaqPjjFhh7G6RwLgJVCHQSCuVTBOIICvUAsVsMUm67SvjzGq7XvtWXOUIgylopjoPL5UQi7LO+ciPGOKBSKwYSdhDwEE37SJhyh2I0M2AD8a3f/EzN7BNNAs3czpfre4+5XmtkVwBXAy44ZyQq+sr/t6kqHD4J6rlolm04bBCmlQSDqtmgdt/q0gWkc1qtjirEBDJI3FlUPmMQ2f+wks/Y0PWBOUVi7C8oardPMrpNl62r7Ac+ybF01TJEJFZixqI4BaPr2JhHa/tiqBrOryu9UJi+tidYLGGnWYCQ3TDXLcBx2XIDNDW63zZ/fa2YfZ/LTv5gHXUauYnIYOa4Aq4u2JhxWOhCUn5U6Qr9oG7SOk1hqjCRTB8CapF+NOkwVvQ0xd0nrqlFFd3aG9hs5ZZmgILEmEaFmYbJ1MkR9i1Jh1xgTGZFQkQm75jFrkfyRqXp0JM6djKLzK1YWMJQjaCCxCO0Fo6IPcEdYgozorvaAmdkTgKcCHwTO3uY+8jmmEsXjBCjQtzXhcMD79hkwX9kPzTNrC2xjvXHM6RhIegjGXpBZ6yZR1zyzBtYLLgYgEKFAcay1tbuyBEDx5ioov6OWKWffmh7N8VUYRRQPU4YJ85iD1hTVxlsn7Joj+n2phJ2qVDBSeWeomrYsGQWE4x4C4ZpBzHuOXRNgZvZw4DeAf+Hu92y/2+7ubkeR7WZ2OXA5wPnnPgxvLMAM8FEgagSiDphcIBu/cauygNSVqWetNf1I85TKMEqOgaZcEklJnxckpUGApA9OgujuoVXX2bA3xqtRBFcKr5qBybI7voIdl1Nl/XWtkbnJqYRdJHFLrIydBpGkSWEXkOwBk2FmCybx9RZ3/8356c9vefCb2TnA7Uf63nmq9QGApz7ly90XD2+6Ni/9VC7YGK8jlPZ3ZiX9asOmSHwwDbpujeLOYRH98Q9jc9dGAxQFk1Zd41ZYPM6bq2qdpX3/kwyVKY+ZpNxZ5SW37GWYVF3/UzRh1zymsgcsiGV8KOdOYdjW50HUEsTsARNgU6rrjcDH3f1nt33pWuBFwJXz/99x/GAFusY29ICvtM/S1DpMmaXG2LAhcEFUZcCqRIRKjEiUfT8KlzbFxauKhB2iMh6JQ5tqkGmRZddao8vWdSicRmVDyQV/YtHKMCUZy4DCrjWqOWChsj/BRjIo3hIlfYARceLcoDwFdiMD9nTghcCHzexD83M/wSS8rjGzy4DPAJccN5IV6Bu7IAI+rjWPaf0+iamD9yvtxVIt7Q1DAOsHjbDrBadxrSCZrSXI/lSX3OWkuKjnwxEYC4p6wERucoBJatoEjo0U0SiZUbPhUAjxTmNOpOitAzSz60S9Wi7KWSqEnWweHEiEnWwTG0jYRdI0uixgPLIHTIC7/3eObvTyrBOJZRhmbf8J3u8HgQDzftBsZvuV9kGHDaxvPwvNV6pkILWpBJjizbAXuSCKMmCS0rNeFFdSsuCYQthBrD44RUnu0El6DDXjCAwX/I2VXrPRUGXrJBt60QZZIeyU2TqVcYqESMIukKaRiMWQOsayBHHvY1jrEsQ64ItHtI0JeB1wRSf7sI41ziq5yoofJKV9ksHZVVMaBGiMSFYEm/lhRFEiZqW9EcmEqAxTMFvLGCUlR5JzVjCvS0mocQSiPktFtk41F1DRWwdZhgm6PhoX2cVH6dkD0bGVXBeD3OjbTpYgBsAKpbQVYJU1TV9ZvzrNq2qNorRxWNdY8ddR8warcIKsVXLH2xQDuVXlkjJGiVKQjA2o8+y21nEpzW+cTHHbY6awB9jqexGZvDTGikvGEXgpU+zWcQVrtQ4UA/wU2TqIVYapFHYSImWqImUBT3/N8ZDJEsQ9j2GNs0qFVUaFUOr2A4eah/V+X/NslfUiK35F/xezbX5jcw+vo2SDLBmeXcq0i2ktGFXvDrUikQoSG2eRPXT15ntZd5vKxBoz9YAp4k5tvE2pBUXWFjfN34NixhzoXFwlGfGSZZgQrr+uOYH6yoAwPXsRZYyTLoh7HoPmPWAUKIIMWPVBY8LRrUFpfBWvHQxt56sBuKiXxldWNUYkkr6q9k6QVscpbmtnC9mdw6rpAVMIZpURiWAenOG4YBdjdcQHgXsptX3fXhlBsFavSMStbCC3QCi5a8xYvFZJyaQG17iXuqjMN1h/nYRIwk7xHhNRx6j+HvYYoQUYZpTGAqwiEHXQvldtxntBueSwBiua5niZFT+NXRtrFVm71+ZCwel1a1Whit06bgVFtk6zkdWUS049ShqHyaXvgxManLTGGDXjCIRz5pojstTzCiYomZTMrlP1qxXNoO9Iwk4z7DxgDxhZghgAw6zthaaAxCzDvNf0lnXtbfgBXNBEYMO6zoq/NXXEBVb81FGzkesFv6+tzJoCyew2az9nrjoS8aHIgJX5Utv6VybLAoIVUR9c44u3dYqJeNMoAkmpnKDawKtoHEGtuCADltk6YRmmgoBjDiTEOb2SUyS4AGuPWafJgAliApOPcWtqPzWCNcb7fVhtvUEeJQLMhg2ZsJNkqxRjAwahxb8qbvOYrrP4b52kcJ96qhpfwL2iMbZAU+dvpQqOQaA+ONCUjJaqKXEFyd+YStRJzJk6wgysUpmm5JiDiebCLmgiKXvA9jiGoATRh+YxYc6qSe5GKRwbkdhouQ9Tf1nruIPAMKR0FEFPkdVRY0aiGJyNyFVPOWOtNaq1MgrKrjS/MauiIpbBUfgrSsoaS6A+OBDNmNP0PzmdxAlSYnDSIRpHoJkzp8jWmahUUFmGqXhPiCbsQpE9YMtJsV5y+ndllWoaUdO+XHAV7wRXrzrgimOwcobAhKPDVX1wkmHUggxYUZUciUobFTFV5TYKVOWShkYoNY8IVCT9KRPL3Qfn1TDJrjOQa2U1WbZOIkIFvXVeLVwZZnNhp5pZJuivM4VpzA6QPWBLiiIDNvooiVutb17eKJuFVjcFjo3zgOvGpZhWO4kVv5VuyjA2RtGvZoMmriqrJMmsVTR3vPuu/eaoIMrWuUQtWY+gX212mGxNRZOpwjWljQozA2EfXHtEm86CxOBE1WsbqgwzmKhTIOuvC0ba0CdNMetwV8xQaf8rlM1CK4v2G66ux/u1xkGZylc6kRV/67ltcx9c67ll3qOZhdYvBBnLecC1Ysi1glqR1CT3is2h0IhEgWgfpxrILSkb7UdRZk3TBydBdAddYXAiEXWgy9Y1j6oRdRBQ2C2B8DgunhmwpCGq0kaJuYdyFprCNGRsL8AcMIEI9To2zwI6UOqIN7biV5RKTqxp6rslWTXN7DrJ4OwtFJk1CaOkVk5RLkmVtMVOgkZgly5xmKwmMuEQ9cEphFLVuVZGKUXVjiNQ9JZFM01p/P4VsgTR8CCmNKdCCrAdRGLuYaI+JdUstNaixgdcNmNN0Qc3aoxI+o3mMWVEcoIUWm83H5y9Fbf1datWzayqvtOIUBdsOIqor2pA5jDZGi+aGxKTuYniGCiMF4SDvgUGFN7canVC9p7YeLmq4eEQpwwzbB4pM2B7G8epjW9LKkSSktZz0EA8C635+0uPCzJgAK4ow+xHBD4k7edfgaSsEYDSacqO+oVE2LlieHjfi8RHEHMTgGFEkl5TmHCorPh7UT+JYq0VXCDEVX0vshImwXo1fXDImpVU4wha/87MNDPmQDdnLpnIHrAlRGVDH4mchSaMWzq8F/SWDe1jGuuiWWhVI+yGzfYXW1kGTFTaqBgqVdAYkRRH0genyNaJjEhAl60SRJUkbV1hxgIiK/7sgwONaHYETpAi10rQGZzECCome8Ai4BJjC8UQYhUSZ0XhLDRrfMp5HSSz0LBhcm1sjNehvRMkaERd6TT3Zeuo6QFTOEHWKrJLFzpBtkZ5VzbKRVZo8R+mBLEisTW3gfZDyQHZ2IDsg9P0wSkyoQVRv7FuHEF74vWAefaABcDblyACmAvK+oKJuiiGIdb1gjloUMc16Pc3j4sP7fvgAO8FlvkDuMAJEoETJIhmodUxjhMkiGIqBZhCMCo2HFXjMDmMGhEqEnWa0mHXZMBEA65tCOSGqeqDU5yztWjGJ6iydc2jon2vDYbkXDC7CHgd0y2fN7j7lYd9/V8BL2Hywf4r4J+4+2fmr43Ah+eXftbdn3uq64mjCo7ANCug7WbWrF/6rBoEMgzxQSPsSo9iYJebog8OUQ/YPolYNFEPmKRXC90sNE3pmeYYSEoQa9H8LUis+G3uWWtMMYFpiqhcsiDpX1XdnzeJENeUdyrKGkFZxqUqbVQMjm4eEitV1gcXI6gYb98DZpNhwi8AzwFuAa4zs2vd/WPbXvanwIXuftDM/hnwM8D3zF875O5PabmmXVMF88G4HrjV3b/TzJ4IXA08FrgBeKG7H8farU6Ziqbr6iVDiM27UFkwBYpZaMV6iWGIbBZav7+9qKmbEidIYw0UpY0CK34QOUHWIsuAuSILqHCpU/XBgSgDJrqLrOgp6hGNDVBkQsnMGqLyTkFMIFQfHMU1WQ+BUPKqy6wlEwITjqcBN7r7TQBmdjVwMfCAAHP39217/QeAF7RexHZ2UxH8KPBx4JHz458Gfs7drzazXwQuA15/zAjeXoCVblVS0uauKeuLJOoilTbKZqEBTmPXxtJLhhTpGsMDWfHXUTivK8osNE0fnG0Nz25JrdAp3g9UTpBosoCKc1a1Nywicw/Fsa2gKRXUuL5JhkaL+uCoYK0dTN1wwa7DVL1l5oLS2Xg9YCDJ3p4L3Lzt8S3ANx3j9ZcB/23b41Uzu57pXftKd//tU13Qruzezew84DuA1wL/yswMeCbwvfNLrgJezfEEGDTfeNZxTbKhrzZorl+BBBgEKm1E1LNWeqC9sHM71DwmveZmBP0+TSmTwrFx2BCZe4hcEAWlgtNQX0HcftH+OJQiyqoBvei9VjITT1GG6VP/T/O4psmsSUq9NQO5Z9uBtiGrScYGTGuNkv1xmcW/RNghMDgJWIJ4CiYcZ80iaYsD7n7gRIOY2QuAC4Fv3fb0V7r7rWb2VcB7zezD7v7Jk1nkFru1e/+PwI8Dj5gfPxa4yx/c7d3CpFaPiXttPgPKZBmwQZP9WXLDEMhZaABeFH1KaAR+GZq7NlodJGv1XlTW1y8kWUDZkGtBbxlFlFmTWOYXGAQ7b9GMIk3pmWggNx7GiGSKqwnbfIOgyiwWTdGoRNTVMmWVBCgEvlObC8aA+utUbOjvcPcLj/K1W4Hztz0+b37ur2FmzwZeAXyruz/gbubut87/v8nM3g88FYglwMzsO4Hb3f0GM3vGSXz/5cDlAOef9wgY2975d6CqZmAp/mC9by5AchZavFloJiiXdFVpow/NLwpeOo2z4iDKrNVRc2EsXXsjDkWpDfPGIJITpEIsedGUYUpKJotmxzCMOlHTGnfJPSlJuSTAEGcgt0TUlREfJDMONP2FRZFZi1mCKOA64ILZb+JW4FIerLoDwMyeCvwScJG7377t+UcDB9193czOAp7OZNBxSuzGLvvpwHPN7NuZarEeyWQLeaaZ9XMW7IjKFGBOJx4AeOpTvtwZGvfT1EGS+ai1n7IfjTFrv0ku1ocqbcxZaMOUuW1s7mGAdwIr/jrgihrEfpwyYQ1RzFcDZFb8mh6wUdgP2Phub62hXCu3YjelE6mZWqfuh9YojEggVrmkaBi1ZJRScY3DpJmuD06QYQ3TBxcyBda+J9LdBzN7KfAuJhv6N7n7R83sNcD17n4t8O+BhwO/PnVGPWA3/7XAL5lZZbpVcuVh7oknxY7vst395cDLAeYM2L9x9+8zs18HnsfkhPgi4B3Hi2VescYliO49DO37adz6yYK8IVZEM7DI0sZohiHmPdbYfKCOa5PBR2u6VWhtRAKi7M8oseKnX2kvwOqo6dViGh3QHIWLWBG5VoImA1ZEGTDVEFNJBgwk6kORpSloRKjKrbDX2LBryJJJp8QZTi9GIcbd/Z3AOw977pXbPn/2Ub7vD4Gvb72evbQjfhlwtZn9FJMX/xuP+x3u2NB4AG0Z8a795tDLonmZmI8DXtqXnlFmIdqYSiwRpqBY33x4uFkfaxaaol8NmveDwtYMLME5W8f2e+TSTdb2jQWIzbFb40VXhqlAMmNNUSqoMiKRuSBu3VRuHVfVAxYns0ZV9SoJEI05YBTM2gNUoi6t6JncK5dAiO7qbtjd3w+8f/78Jiaf/hMIMGIbB9suqnSizUYvedOqArFortnQd9Y1Fx9KQScpbURzLVDNQlPcPJU4K4Jkfp/bMJULKmhswmF1nESCIgsowEAyC01T2qmZheaIHDEVjo21Ikn/9L1u1ENratU4TA4jYSaMVdFKRWWYkpllSofJ1gJfZEKixNFkwPYasdMRXrHNtgLMSkcVbGa99JJSJsUAXhdlVIqlYUjE0sbWWFlF4UPvjQ15Hogr6ANr3asGOtMUlQBT9VVN9vaC8iDFOAJmp7bWcWWCRlaD2B6VEYniHbzvdOMTWtOBahZa67uTuvJLZblkY4LqmJO0oQ9FnJ3rETB3ykZjF8R+BesalzUCNqxLMmCSsiuRFX/1nIUGOQvNSq+5MHb7oW62jVmAqrC339e+VFBlxV/QlPX1i/YOi7ViIidIGxqfWyArl5T1qyms+GvRlAVIjgEivagqw1RkgwljxT+5CqrQOEwmEzUzYHsc9/YXxfkC3hrr97V3aCu95K6/1IpfkLFbdsMQmGzz3dtebIv1wWah9e37teqgyVx7r7Hib90Ty7wpUGzkStd+C1M6oQmHwjRElwWUOGKqetZavyGoxgZArMxazm2LUtgpjxyK7AELQB2xtfvbxuwXFMmdsw7vGg+fBUaFTbbQil+BrFcroAhrHzPQLDRV9kdw0wDABU6QJng/kPRpbdE8CzhqxIcoU6Xqg7M6StwVVWWY7X9fIvHlVdNfN4gGgSks/guacsmBubyxNYr5YppqbwkBdUz2gEXAHVtrmwGyftD06JRO0sxfNhoLUMB7jQBzUf9T9fYb5AKhShtzFtogGUY93YwQXGl9YBqD2DhsL+oBU1jxe6cbyN1YMJkgJuicICV9cKUDSRmmKlOlmdMkK5eUDGKmfdwqOgh91cyDEyWVZAO5EyAF2N5HUILogA0CA4qNNcmdTl8cbN7vMA3gFVjxC2IC1LI2Waa3jEmWNioNQ1r3GFrRDDqnaBwmqcIyzMaEsuJHM7PMQVJ6ZooyTMBKwQUZMMmcudJeiAM6YSepQMxyybnsRBNXgiazlkykANvr1AprbTf1priDDNPdQ0WpycrB5hsOr2NzQQPg/arGiKQf8LH9bK0sbdQQygmS2bWxMV4GzcZAVYYpKnXWOIkJSnFB1gcnQbLxRrPeUtuLRUROkCojEpXjm0TY5dw2IE67VkQd40ZNF8Q9jldYb5xVGds3x29hjevnvRRKv9I0JoCvaHpeaunaZ8FKT1UYe3SrjIKeta5oZmupRN2yz0KbBpkK4narmvVKbOg1RiT0+zT9ZZGs+KH5PLgppkIsajZEUw+YoDpE0a+m6K0DQDQLLZIbpsqKXxVXgQEeRdnpcEgTjj1PdVjbaBtTNHTVSsEbvwkYUPpF87uHDwx0bYxt7JMcX188onlM0DlBphV/HMMQ1Sw0N82NHhcdA8W8Qaqm/ymUFb+q/E5wY27qrYvjBClZK6K3WoWggcysyeO2RuQwGZAsQdzrOLDZ+I+guuiu0RrW+k22FNhYb9/IvlIpgruH3q+0v5Neenzz3uY9Kl43NUOuXWPy0iEoERMRyTAk3Cy00v7vdtq/KISdYl7XOJVLNn6fCWXFD5p+NSbTkOZxS5H1wSnQrFWUoREM+QZ0PWuRessUoq6O7cMG1TEpwPY6Dt64os16Bxpn1bZo3V9WCqawxK2j5OJVSoevCDJgK2e0j7l4OHW4r3lcs54iEHZVUHYVqVetWM/YeA7aVtxws9AUCM4v71bbxy39NG9R0BcrQ2VD3xhHY3BCqZrsYpSyM5hdEAPNQovUn6Ow4gfizG0LiOcg5r1PBT/U+ELbu6R51+oG9IrhhYLm5Y0NiuritdFWMfs8X625E+SwTlXMQhM4AMLUW9aa0Ue6EiezpihrnOJqSlFbx3UfJIYh2KDJ0lSNEYnKCVJRlq2w4oetaoPGIkxlGAKanjXJzcnS3HkZRIOzVeWttaLZOg4awTgMgaz4aS8WLZ6QcSwzYHseN3yz8T9hcKbhGa2p2CjYxhxsOwcNgJWFJgPWL/C+7cXLSoevrDbfHPm4rhlqK8oqjVUw1Nf6UFb8kUobvfR46zI56zXupbISxL593Drg/b4pC9aQKfsjEAkiK35KJ5qFlk6QUiHaGlmpoIAsl0Rm8Z/sSUILMK9GXW/b82Cd6uQfYBBk1jqF45emD876Hlr3lpWCrT5M0Ac34ouDjaOCL9Yk2YTaP7x9pqYAgXrLFEgt8xUzsESz0BAZkTSnALVvL2q2ShsbM904iiGUpllomhJESRZQ0a/GZBqiiKtANpJBQZZLIjEMCZpIygzYXseNcb2t25NZpYyC7I8bLhB3pbQvh7DNAWp78QFgi/buXKVfNH/z9mFDcvHyfh/eN29cpK60z4CVbpVRtEmO1F+mINosNIlY6vcL7O0Xc0xNZq01Cgt2L72kZ2+ahaYQYKNoaLRmIHfrcTIAPmzGcoJsHhXh4GyFqVr7kA+SGTDIHrA9j9fCeKixACuOj4I37qHDRBb3rbF+FPXB3Q8r683jFkGtv6+IMj+loy7am4bUxcObxxz7/ZINvUp8hCptJGeh4YJMFWjKJQGXZBc1Vvze7WtuxOGdcLsg6ldrbkZSOnxob9JlpYhEqGYOmGLUnvXEmYWmegOXEFPIZAZsj+NujBuNSxCLSwbAeTWsdb9acbD2BW3WjxSJE2TFNhtn7AqwInArHAaZFb+ttD22Xjrqxl1NY06BB0aBqYOVXmIakrPQchbaVtzmzKWN7eNqDE6sjO3FeB1A0BMqM4soo+jYipwgFQTqr4uXWRPEVCw1oI5xTwG291EIMNMIMNwkYklB6TSZusIGdALBuGhffscgsuLv+/Z3UEtH3X9P25iA10Ey5LqOGsOQaLPQVKMDWsdVODZCzkKDKaum2XSKNt4Kt0KRUNIJO0GFTKB+NcgZazIkmbWIQsayBFGFmZ0JvAH4Oqa/uX8C/AXwNuAJwKeBS9z9zmPF8WoMa417iopTBkHJ0eaAlfZvLz52zd9gRnO61tk6oNvs25uclErH+pQNbIitrEvKIQqC8sZS6FZWm28M6uZ+RsXMMmGpoCKrFKm0cfRR5trYGtUsNBQD1FXz1QA3hZGSakMfx4pfNuhbYcIhGBvgMhv6MZ0glQTSi0oyA6bjdcDvuvvzzGwFOAP4CeA97n6lmV0BXAG87FhB3I2NjbYCrJhTBGYZdbNrLhIAeoVhSFepEhHaY4Lsmg9dc3Fri4GyIaj1r3dhK61LUQudoGetrJ4hKesbS8/YOmtZpgHXrTM1xfpQpY05C22QlCBO4kOQWVM004AmS4PGit8BTGPx3xpDYggKgyZPIcusCY6B4j4EiG4alE4yD05CUB2TAkyAmT0K+BbgxQDuvgFsmNnFwDPml10FvJ+HIMAGQQ9YGdsLsK4T2Lqbg6JcsmgyYCojEkVmsSyGqWy0ddyyjq00vtIUKPu+2DYm4GesT5m11nFXzqAKLL3r4szmMb1bXfqsWs5CE85Cq5osoGoYtSrzEaVXS9enFKdfzRGYmwBeShgnyClulMxaPCHjpAuiiicCfwX8spl9A3AD8KPA2e5+2/yazwFnHy+Qu7G52VaAleKYtf+DHUv7LA1AGfrmJYjFnCqYWVbHoinDFIjQbjFQBSK0HzpsX/s7Z93iruYxbX2NTmBEUlf2w9DYDbP01DPam7G4D0tvGFKsZ3TBjZOchTZl1STzCBRrRTJfLNoMrGXvV5viakSNRoQqXBBrnBlrEXVMmnBIf+Y3Aj/s7h80s9cxlRs+gLu72ZEdK8zscuBygHNWV1lvXIJo5hzlR58SRRVX8EZYilMEF9rS7Wt+DMycofEsOIBSRrrWmSpgcf9+St8+7r7Nu5rHLPvXKGufbR53MWxSH/7I5nHXV9vH9JUzGRQCrIeuxDENkWQBEWmPSLPQVL1lveCcrRrHRln5nciK32g/SoXSiYSdZmi0rARR0rOnKUFEUIKoybBGFDJpwqHiFuAWd//g/PjtTALs82Z2jrvfZmbnALcf6Zvd/QBwAOBJj3iUt86AARKhpIgJ0KkcCxXCTiRC+3WBWCwjnaBccnFwtXmPoVllXGu/me9W11m5p31WaXH/zXQPa90HZ6ysPqz5Xcm6/0yGfn/TmCCaq4VO1OUsNN0sNBrP13IfRJb5pBW/0Iq/tVCwOkoylqoyTEW/2gNlgoq5poGcIKPhgCukqNlFTB4UHfAGd7/ysK/vA94M/G3gC8D3uPun56+9HLgMGIEfcfd3nep6dlyAufvnzOxmM/tf3P0vgGcBH5s/XgRcOf//HceLVTE2RsWg2Dh/AiYwyzBzyTFQCdtuUJRHOZ1AhHbdSBGUYa6vtb/jvVjZZHFP+wHP+++/i25/W4MTK5UFf9q+Lv9Rj2JjbH/He/Psp+Krj2seV2HGoiRnofXNjTjMeo1hSFrxT4is+BUuiIpjoCrDtNJpji1o+rWiOEFazExS6xJEmy42vwA8hykRdJ2ZXevuH9v2ssuAO939b5jZpcBPA99jZk8CLgX+FvB44PfN7G+6n1pt/nHfTc3sSYctEDN7hru//xR+7g8Db5kdEG8CfoDp/to1ZnYZ8BngkuMFcTfWg2TAIqEUoKpjG0UwFnNMcPlSiMVSKr2gXHLl82c1j2vmPOrW9kYkq4+5h9Xb3tc8bvf1d7D52K9oHvf+R32Fpr9MsqEXuYgJ5sxNmY/2fwvWuhcSsDqI4o6iY7AhMXWw1vMWAWoVxR0xQZZGUSZntWqE0qCxQZS4FdYqyaoxCo7BKBK1YgQliE8DbnT3mwDM7GrgYqbkzxYXA6+eP3878PNmZvPzV7v7OvApM7txjvdHp7Kgh3KlvsbMfhX4GWB1/v+FwN892R/q7h+aYxzOs04sjrHeOAO0DHWnx0OR+t1CcXxVzZpVcBzcNcdg9PZ3+MZqbAriHhp6xsbHoLrxubVF83Phrk3jE/e2v4D/ebmRe8ZPNI3pPnJo+MOmMR+MLSjhMU0JTzGBWERThqlaay8oRTUKReQIWmgft0Nwc5aezh7WPG6hYIJj0AuyttNa2//tRlproVAE1/FOcA78+drvNI+pxxT7unOBm7c9vgX4pqO9xt0HM7sbeOz8/AcO+95zT3VBD+WM/yamNNwfAo8A3gI8/VR/cAuqG2tj+zdZD5IA29rIKwSIQnxUN8mxVYiPCoy1fdxNb//GUt1Yq6V53M1q3D+0PwZf3CisNb4p58CN97fPfNxjB7nF/rJ53LvXPssw3ts8rrsg+yNDc+NEsZHFCqbY0FsvEaIbklloHUVQglisl/zOOsnsuo5OIJpjiVCNEO9E54FC2IHm/JqEXdv3gypw9VZzCjeqzzKz67c9PjD7RuxJHsoZtAkcAvYzZcA+5YrboSeBY2yM7f9gNVma5iEBjfiY4moEmOLEGWpp/jtzYFMgwNbGwig4F+4dOmrjuBvVuEtQbfNX6yMHG5daOM6N5S/xxmfYut/HPWu3NI0JMNZ7OcXy8dMEhQ+94j2xSOJqsoCabIqZJptgdBTJsVX0P6kyKoHWap3sPGjNlpgpoixY85hemseNWtN1kjeU73D3I1XXAdwKnL/t8Xnzc0d6zS02NRM/ismM46F87wnzUATYdUyGGH8HOAv4RTP7x+7+3af6w08Vd5oLsIoooyQSHwqh5G6SuIMo7sYoKL9zWBcIsEOjsdl4bll1uFMglNZG546N9oFvtzu5v2vvrnjX5s3Hf9EJMoxrjPX+5nEn8RVJgCku4wUTxW2Nodl0QpGIMFVMRdwiE3aaTI3k2FIkIkzy76e9SACaZ36A+aiKht6L1tv+byGqBGvOdcAFZvZEJvF0KfC9h73mWiYTwD8Cnge8dx6LdS3wa2b2s0wmHBcAf3yqC3ooZ+Zl7r6V0rsNuNjMXniqP7gFI8bBxgLMXTN/YBRmf1pTXSPs1qvCfgIOCsrkhgprggHP92/CeutUFfBXG+37lNZ84I7S3tjiTm5jYzzYNKZT2RjubhoToPqAu8IsIpL4Ao2oMU1WSVWCqBI1ko23ovyuSMrvjCLJVknWaiKhJBL4qv6nzgVlmCphFyRTBTpxG5HWbTBzT9dLgXcx2dC/yd0/amavAa5392uBNwK/OptsfJFJpDG/7homw44B+KFTdUCEhyDAtomv7c/96qn+4BaoMmAK8TFWlQDTZKoUcTfdJOV39w+CtVY4KDAkum+obDZ2T6rA3X6I2rr8zta5l/YCbG28m816qGlMp1Jr+2yds4lmXHAkRHdQrRClBHHayIoya1EyYCqRoCppU2Tr6EL1VWlEQicTShLBGCZTpTH3iGjq7Wgq0dz9ncA7D3vulds+XwOOWN3n7q8FXttyPbsxiLkZjrHROAOkyv6MrjG22Bg12bpNwR/t+thegDlTVqk1m+4cHNofhHuGTTYbb+gd557SPvuzbodYq/c0j7tR72esbW2y3esklhq3p/o8EjIOWSqo2sxHKuvTZNU0mSqd+NCsVUHEvqrWyEobVb8ziVg0QdyIJYiaSrS9RmgBVh3WGguQiiZLs1kltpqs1/YGHxVjQ3DTf22EsXHcCtw/1OYGFJvu3D+2V3b3cIhNaxu3UrmXL1Ibl7UNdZ1Dw51NYwIM4/2abJWPZLZKVL6iEglBsj/RjC0kphaiu/6dLWR9Va0paFwQs69K4wAIujJMhQ19kd3siodqvNBeIrQAc59c5VpSoblJAkw9RYoMWGs7b5iydYq4a6MzNlaL1eFewfDCTa/cR9syOYB7yt1sWtvsT6VycBQIpbrO5tjegKLWjbm0ryFeiddXpbJhj5GtQmW+kBbsslJBlanFsvdVKcoaixd6lQFFoL6qOJmq6b27NH7/jihjnOW4lRpagFWM9cZiabIf1/QUKTJrOgHWfrFrY2VT4Me/5kPz/qdNRu4v7cXHmt3PZuNZTU5ls7Y1tQAY6wbVNX1VacEeSHyJjCJM5QAYKFOVFuxaoRSlrypSnxLEWm+kTFURCLCQeGbA9jzONKuoJdVB0PojFGDtg04CrP39h7Va2RRsvO+39sNnN9lgzQQCzO9laNz/BJNdemuqb0ocALNUcIsYBhTaGVgxShDTgn2OG8yCPUxflapMLlD/U6QeMEWmCrYEWPaAgWYe714jtgATlMqpBNj6SPPyO4BDAgE2uLNW2wulg77BpqBM7KDd1z4DZusc8vYGFOvjPYyNRY17ZWzsKjjFHcAFVpDTtD1BXBVpbLH0M7BUPWDEsmBXlMpJyu+sk/QpGV2w7E+ctfYuGB+BpipAlalSZNaiyhgPu/KHTmgBVqG5WcRYkZTJ6QTY2HwrO3hlrXWPDlOmaqTthr5a5aDd2zQmwKavsV7bDwverIeotb2ocR/w5uK2CmJCLPGlIWdgQRGJBFVpYyfoLVMcg2Iaa/tCl5kqUV9VrOxPpLVGylTp1hsNTxfEvY/7JGyaxsSbxwTYqFUiwDZ8bF7MNVJZFwiwdYUAo7Lh7fufBl+XlArWOlC9fdzmphbQ3NI9JoauVDBGCaJsBpawVLD15ltlwa4SHyoLdlm55JJnf9IBUGUeRJhM1RRX1bUXj9bO1nuR2AIMGBr/lgafxFJr1mulfa4KDgk23iOV9cZOfTAJsAFB/1NdpzbuLau+ySAo66u+ThWU9elMLSK9CwYytghUKiidgaVwlBPZhKvER7S+quYxReWdkfqq0gFQI5SmuHEyVYX215uoeaQsQdzjVHfWmgsw1wgwH2mdq6pUDgmE0qZtcMjaZ5XW/T42aW8WsT62LxUcfb35sGAQWbBPkQUxI4kvyBlYIgt2K5pyQaH4aB4zmAV7pL6qztuvFQg1q0q1VlVfVScSYNHEUms6hYmSxRMyTppw7HmcKbPUktGddUE2YZ1NiQBbFzkAKsr6NvxgewMKKqOv443L5apvyizYNaV90cRSa5TZr0gzsESZtSA9YFulgq2PbbGOIhC3HYswmSqVEYksqxQsU5VlfZqyvj5QZUT2f8341GJ0uhNegG023syOXtls3KcEsG6bzfufYBJL1RofAzYZBH1Kg69TBbbmikzV5AAoEEpe8RRLoqg5AyvSDCyVBbuqBHHZhwWrsj/pAJhlfaA1oJCIJVFmqXXUqJKuhl35Qye2ABNkq0ZGiQHFJhsM1laAOZUNa9+ntOlrGgFW2wswp05iqTHVB0mpoC/NjPfjEcOAImdgqUsQNeK2ecxg87okJYii7E8sVz1NWZ/S1EFBJAMKSVkfmmNr6IRdJJwcxCzDzP4l8BKm4/xh4AeAc4CrgccCNwAvdD92DVjFm4ulkVFiQHHIDjJa+w39hh+kNrYKrz6wWQXOgvWQJANWq6JUcBQZW+QMrJyBNfVqScq5ZD1gmhJEifgQiUXJDCxbhHIAVPQq9a7ZhijWChphp+ipgliZqi5QuSRAJxBK6YC4RdrQSzCzc4EfAZ7k7ofM7BrgUuDbgZ9z96vN7BeBy4DXHyuWm7NubTffI4Omr8rWm/c/wVZZn8IBUGBA4ZuSbNU0q6pxVslr+5hTYEHMWESagWXWa4wtRH1VqrJG2byqMOV3GhOOdADUlvU1j5mZKlmmKkpJ3xaKY2CCY3D6y5i47FYJYg/sN7NN4AzgNuCZwPfOX78KeDXHE2B4876qwYbmpYIwCaVRUNI21PZicfRBZJWuiYtE1DkplpQzsFqjKxVUzqtqHlfhqGd9KAv2WDOwumC9Sstegqhy1ItTJqewSp/iijJggUr6MgP2IMuw+9pxAebut5rZfwA+CxwCfo+p5PAufzA9cgtw7nFjCVwABxskDoCbrEnEh6JXq/qgsWD3ARdkATWmFtH6tOLcQUV111+QqVKuVSE+OkkJYhwL9i2hpMpWNY8p6quS2JqnA+B8dmnWGkosBcrUqESNpAQxkFhUkjb0Iszs0cDFwBOBu4BfBy46ge+/HLgcoLNV1hubUDiVDQQCTNT/NAr6n2Slgr4pyVZlqeAWMbJVqhlYKlOLUjR9VRphpykVjGTBrugBKxTJvKoeTXYxkrBTmFoYRi/aeneSktE4BhSyDl7THANFWV8hjgiNKmOi3QI/GXajBPHZwKfc/a8AzOw3gacDZ5pZP2fBzgNuPdI3u/sB4ADASvcIb91XNbIpyVRV35T0gClEnfsgEWD4IMpWRRRLrdGU20SbgdVa1Ojc7zpJtko1/ymKBXuxTlR6poobRygpM1VRyvrCZaoCOfUVdGKpeUxTiVBB0KCkC6KGzwLfbGZnMJUgPgu4Hngf8DwmJ8QXAe84XiDHGWhbKlcZJUJpdI2w02Sq6mxs0ThuWrDPxCi3iTYDS2KVLurV6qxfegv2zhaauIJMldKAQtUD1jymcFaVqqyvNapMlWFhMiqZqdoSoc3DinKL8YSMe5YgSnD3D5rZ24E/AQbgT5kyWv8VuNrMfmp+7o3HjYU374GqPkocAMe6obFgl2SqapYKykhjCxP1P2lMLYrG1lzgAFisC2XBHsl8QWnBHskBUGGXPnftSeK2j6nZJPeqWVXBMlVxRKiqt655SKK2lUXb2Z0Mu+KC6O6vAl512NM3AU87sTi1uQugUxklxhbCvqrm1CwVBJZ+BpYtZHOlophwqCzYVdkfRbmgSoR2LCTZH1X5nWJelWpWla5PSTH/qRPF1ew8FWJJVdanMIoATVmfaq2yYyAIO4nb9nEjkhmwPY8zNhY1zkitmlJB97YZIEc0q0o2AyuJViqoGm4sMbaQiDpVD1j70sZCJ8lUqfqqlK56zWO6blZVnOyPsgSxPZHmSqnK+kDXA6VAImpE+3hNGWb79UaUMcvSrBJagDnevKzPqbpSwcYCDJDEzBlYELFUsPUlXNX/VIqo7EpyDDQW7KoMmMqEQyFue8Hlx4Q9VbohxDHMF3QliGlAoSrr60OJUFUPlE4stSZNOLawNOHY+3jzHihtqaBKLC3DvYJjkTOwVDOwNLbmfRgL9k5UhtkhMqCQiJpOUn4XzQFQYZeu6qladqt0iJWpUpXJRRJKKvHRqTJgIhOOFGETy7CrjS3AXJABU1mwUzUZsDS2ECEyoBDNwDJZ/5MqbuNMlcosQ2lAocisiRwAVQYUCgGmEErKTFWcEsRYc6VU2R/VAN4otuYqkaDqf1IJGpUIbf0XFlHPOTtvQ29mjwHeBjwB+DRwibvfedhrngK8HngkMAKvdfe3zV/7FeBbgbvnl7/Y3T90rJ8ZWoBNv6TGPWAiC3Y8jS105AwsSQ+YqgTRNJk1lVW6alhwmP4nafldjGOgGuyrE3btUc3AUpX1KUr6IG3NFSIBdOYTsrLGQCI0InXnt7ZXAO9x9yvN7Ir58csOe81B4Pvd/RNm9njgBjN7l7vfNX/9x9z97Q/1B4YWYODtjS18IEsFlcSY96KdgaVxFmyNslRQJZZao8qsFVUpqsgwJIoDYHELlalSlCBClvVBvAG8Ucr6ImWqlCV9UUw4QqbA2JXUwsXAM+bPrwLez2ECzN3/ctvn/9PMbge+DLjrZH7gaSDA2rsg6koFlz1bFcvYQiOU2s/AUpUKFlvQKfqfVAJMVCqoEAqdazJgGlETywFQlalKAwpNWZ+pBFigsj6dBbskbCixFEYoodFKEfXXLg1iPtvdb5s//xxw9rFebGZPA1aAT257+rVm9krgPcAV7seeaRVfgNHYsTBLBYXEmoGlySrp4raPqTPhUJX1tUbpAKjYdqrMMjQGFJr5T5HK+lQ9VTJTh0ADeFXiQ9ID1jziHDePgSxbF0XcR+Uk0yBnmdn12x4fcPcDWw/M7PeBxx3h+16x/YG7u5kdddNuZucAvwq8yB8sw3s5k3BbAQ4wZc9ec6zFBhdgBDK2iEaQUkE0xhayGVjCvqr2MRcSo4hCJ+qD02Q+IjkAqo6Bbv6TYgBvrExVjFtS8TJVy25pDrpjoCJOaWOcDFhUTtKE4w53v/DoMf3ZR/uamX3ezM5x99tmgXX7UV73SOC/Aq9w9w9si72VPVs3s18G/s3xFhtegGW2SkGWCpow+xMlAzY59UWyYBf8vkTDgqM5ACpETU8XJlNlGL3EFTWOAUW0TFW0/qfsAWsfE5RW/Jq4rYko6nbJMeFa4EXAlfP/33H4C8xsBfgt4M2Hm21sE28GfBfwkeP9wOACLI0tcgaWyILdCqVo+qoU2SqJUYT1ElGj6gFTZKpUFuwdGgGmmv/UK0pGJUcgVllfpAG8Sqc+TV9V+5iQPWAQRyiohNJW7Agxk4fMlcA1ZnYZ8BngEgAzuxD4QXd/yfzctwCPNbMXz9+3ZTf/FjP7MqZT7kPADx7vBwYXYAqWPfsF8WZgadYqKecSZtaax6TTCDBE7ooCW3NVCWKhUARNxpqeKpWzYCwDikhlfZGc+vpAIgFi9T9lVkl3DKKI0Kj4Dm/F3f0LwLOO8Pz1wEvmz/8z8J+P8v3PPNGfeRoIsGUXTMs9A8tEdvGl9CKhpMn+SEwt0PSAGR3Fl9sBsHdN+Z0iUzXNwNL0hGpK2uIYUKiGBUfKVEXrf4rVA6YhivgoIrdCiCMWo4q6GnblD53TQIBFIoixhUrUWI/i8qUpFeyX3oK9s4VE1Oj6nyI5ABZRpkY1LFchPhQr1WWqIpX16colm4fN/ifiZX8iiaVlF6ERcXZlEPOOkwJsx0hjC5MN9tUIpSgW7FvDghXldwvf1zQmaPqfimt6tYqbpK9KZUChylQpDCg6UfmdRHyInPo6080+imKUkP1PurWmUIq31tK4qkuVXVWz0yWIu0EKsB1DVuXdPKJqBpbRSbJVGgOKhSRux0JULhhnWLAqA6YY7DuV32nKfCUW7GlAIRVLrZk2XO3pRJebZe9/0hmRaFAZRWS2ThCU9uIrLpYliMtLHGdB1QwsRQasFM3wWZUFu8Qunk4j7FwgQtMBcLb1UFiwq4wt4hhQqDJVimyCqv9JlalSlfSp+p+iiJpo4iNFqKoUVSOUomarmuOZAVtSNNsNs04UV1MqKJuBpcj+iATYgtXmcXvbJ9rQa4SSoq+qd40BuSJTVSgsFGLRND1gqgyFIq6u9CxOSVs0Z0HFoVWV30UqQYwm7CKtteBhynEh+8BgeQZMyQSYmb0J+E7gdnf/uvm5xwBvA54AfBq4xN3vnAeXvQ74duAgk6/+nzzEn9R66SL3pDgzsKZZVTEMKMw62bBgWQliEAfA3ntRBqzTiA9RBkxR1tdJMhS6UkFV+V2YQcyiLI3OsbF9TKVIiLKhjzZXKspxhRRKEOcc2AnShOPU+BXg54E3b3vuCuA97n6lmV0xP34Z8G3ABfPHNwGvn/+/C8SbgdVa2JmJyu8COQBOM7DSAVDhADgdAU1ZX2tUvVpFZEBR5tjN4wYRH4V4/U9RNp3K0rMoQ22z/E5bfrfsxxZ0PZHRWAL9pRNg7v4HZvaEw56+GHjG/PlVwPuZBNjFwJvd3YEPmNmZZnaOu9927J/SvlxQOQNLMixXVCooyf5EEmDWSfqqdP1PKgdATVmfIvMRyYAilgNgLPGR/U/pqgeZrQONWFL1KcU7thrMWv/O4kmZyYY+au7uobPTPWBnbxNVnwPOnj8/F7h52+tumZ87jgATXGyFM7A05YJ982yV0dFZe/vxYr1sXlV797tO0lelKutTDPYtmGSw70IlwAKV9fXBnPqWvf9JeXc+kvhY9mMAuqxSlEwopFAChVAS3eQRxNwJ0oRDiLu7ncQZbGaXA5fPjwRlfboZWJps1YIicEFU9T8pBJhiVpXSAVBR1reQ9D8VSQZMZUAhGZQrGsCrcuqLJD6UWbUoG3pVpqoX7TqXPaukEkoQ55wF1bHVEEUogSi7GFSBpQlHez6/VVpoZucAt8/P3wqcv+11583PfQnufgA4AGDWefMSRNEMrCKK24nKGmUCTCFCvX1P0dT/FMMB0ERW6R1Fk1US2dxIZkqhuX6pbM1Vd9EjlSBGctWLlqGIIpRAV34X7XemQFOOq0l5RBJLRXQMouFkBkzBtcCLgCvn/79j2/MvNbOrmcw37j5+/xeAYY039LIZWIJSwSnugmKNTTgoslJBSQaMleYxi6syYJ3ELl2SqcLoJTcNNP1PkgwYGvHRm0p8xOl/UpllRCqTiyRqZMdA1KOiOr+i/L4gs0qg7FkTCPzmEcMmwDIDdiqY2VuZDDfOMrNbgFcxCa9rzOwy4DPAJfPL38lkQX8jkw39DzzEnyIwtih0RdH/tJBkq3rBWjsW9IIesIXvk/RVLVwgwCgsBBmwlXkEb2sWIgOKhcqAQtRX1Rql+Fj2/qdomaplzypln5JWKKWoiSNqVKjOr3B42tCfEu7+/KN86VlHeK0DP3RyP6ntW4xuBpamB6xjQWltQy/KgC1YEVmwa/qfVJmq1lENC2VAMQmw5mFlfQmKO+kyF0SFgSuxMiqRBvDK1hpILEXrU4pUKhdJKEE00bwECmGXcCJ6N544u2bC0YKpj6KtUCilF7kgCgf7CkoQNf1PKgt20QBeQdyOONkfE2WqVO53qsxHpP6nSFkalVFEJLGY5Xe6DXIXSChBrJ6iWKI5zlY+yjUhaUNoAYYZpbT9JxTrJWV9fdmnySoJSgU7X0j6qlQOgJpM1ZRbbE1vRZZVas2UAWseNuAMrOZhZaJGtdZIokaZVWpNJMfGLL+bSKEUSyjJsoCasMlMliDucQxrPq+qs17T/8SqJAO28H0SB0BFX9WKLyRmEQuJWYZJ+qp6s+ZufQVYCK6IKlvzaP1PkQSYIqsUrwRRlPkIdAzSfGFL3MbpKVr2UrloQilKKWrUDFi6IO55rLmxhaxU0BYSVz1Fpqp33awqRQmiyoBCM/+pfVmfoSuTy/6nWFmlUMYWoo1hNFOHKJbeKZQmIgkllahRsOxC6YG4rX9ngc6BLZx0QdzzKHrAVFbpPfso3tgwRNSr1dGLskpd8wG8BZbegGLKBLeNCbr+p0giIVqZnE4kZPldFKEEOkOHSKVykYQSiEob24eUEelvAdKxUE2WIO5xzKx5v1bHghXb3zQmwIqvaDJgrjD26CQCbEFpnv3ZylQpyvoU4mNRsv+pC5T9idSrBdArMgnBesCiCaVl7ymKJJSyTG5i2UWN6rhCliBusQT6K7gAo7BgtWnMzhYsXDEDSyTA6AVZpSIxoFiYQoDBorQ+Arrsz0LU/7QQXBFVAixSVklmwhEsoxRJKEUzX1h2UQOxsj/RRE1rZCW+S35cQXUM4kkZJzNgAbDmhhkdul4tTV+VwC4eoxcZUChEjWoAbywHwPYxlUYRUbJK0crvdA6A7VGJpUhCaYotiJlCKdyGftmzP9HOAwWRxKIUTxOOPc/k1tdWgPX0rAgcAHvvRA6AAgMKMxaS7E/7TBXo3O8UzoLRyu9Ua9WcB3FmKmVWSVfOFckFMZJQgnjlnZK4wTJArYkmlCKJmkhmNGrShGOPY1hzF8AVX2GfoARxHyILdlGmalE0cVtfvArQC66IxTRlfTLxEc0BULDxjCaUojjVRXOpUwnxzP7EEjVRBA3keaBE1Q+pIKpYak2WIAagUNjnbXvAeu/YJ3BBXKGTGFAsSntj92IaC/ZFaX+hNSxU/5MqA6bIAkYSSqDq1dIJpSjmCzoRmkIp0gY5WuYnf2c6Ipm8qIgk7CKy00fXzB4DvA14AvBp4BJ3v/MIrxuBD88PP+vuz52ffyJwNfBY4Abghe6+cayfGVqAmRsrjV0AO4okU9Vb0Zg6iKzSNdkf0WDfQOYLOgGW5XeK0rNIQgmUPWBxHOVkJYiSqJn9SaE0xw30O0uhFOhvLNIvaxu7kAG7AniPu19pZlfMj192hNcdcvenHOH5nwZ+zt2vNrNfBC4DXn+sHxhbgDEN921JRycp61uIDCgUfUo6C3bNZiNaBqw1BQ8llFRZJd18MYGwCyaUsvcnllBKUaMhMz9aoRhG1BDnnI10bm1nF0w4LgaeMX9+FfB+jizAvgSb5iA9E/jebd//ak5nAVYwVhuXC/ZWWBH0P+0rJrEfXwnU/9SLsj+KtRqarJJKKCkyi6ARH5F6ipQuiAoiZX+ypC3WWiGzNBDsGARaaxRBA7GOa0ScXTHhONvdb5s//xxw9lFet2pm1wMDcKW7/zZT2eFd7j7Mr7kFOPd4PzC0ADOz5tmqXuUAaCbZfKuyP4oNvSwDFmj4bCfKJkQTSlF6ijL7M48jCLTeSGtVES1Tk6Im1vkVSYBEOq7Jg9STS4GdNYujLQ64+4GtB2b2+8DjjvB9r9j+wN3djn6Sf6W732pmXwW818w+DNx9MouNLcCguVvf5ADY/jLTF40DoCr7o3IAlIiaQFkllaV3V1IoRSq7KmRJW6S1QjxR05ql76VBu5nP4xCHUDcMAg5ihpM24bjD3S88akz3Zx/ta2b2eTM7x91vM7NzgNuPEuPW+f83mdn7gacCvwGcaWb9nAU7D7j1eIsNL8D2Nd7FdGaSsr5F0fSoxBJgLsqAtY8JmqySQiiBrrdM5SgXxYIddEIpRU2KmhQ1sW6cqIgmaEIJkEBrTSbcd8WE41rgRcCV8//fcfgLzOzRwEF3Xzezs4CnAz8zZ8zeBzyPyQnxiN9/OKEFWDFjtfHOcxJgTUMCuvI7xVoLLpkr1YtKmTQuiB5rWHAgkaDs/YkialLQTETZyEUTtpE2nSk+dEQ6DyDeudAaxe8rys2oPcCVwDVmdhnwGeASADO7EPhBd38J8LXAL5lZZTpdr3T3j83f/zLgajP7KeBPgTce7weGFmBGexOKTmRAocuAaf5gFXFVlt6KtYLS1jzGhl7V+xNp0zmVCqaoUZCiRkc0ga8g0u9LRSRBk78v3WzEeDi+w8fC3b8APOsIz18PvGT+/A+Brz/K998EPO1EfqZMgJnZm4DvBG5396+bn/v3wP8GbACfBH7A3e+av/ZyJt/8EfgRd3/X8X9G+wxQEcSEKfMhMeEQZT40FuyaEsRoQilK708ksQixhFK00rMUNfHOWwX5+9IR6dgqSPEhcscNmAJzdqUEccdRZsB+Bfh54M3bnns38HJ3H8zsp4GXAy8zsycBlwJ/C3g88Ptm9jfdfTzWDyjAvrZjwOYMmKZMTiEUVBkwlVCK4gAImn4t1aZAZj8eqCRCstYUNECszWw0casg0u8LYh1bBdHEh2qMxrIT6XqrZhds6HccmQBz9z8wsycc9tzvbXv4AaaGNZgGoF3t7uvAp8zsRqZU3h8d62eYwWrX9oTtTGNr3hVNT5GuVLB93FizqlT9T3HKJSHWhj7SWlWk+Ij1+4p1XOOsFVIkQKzsqopo5217Yv77fRcmMe80u9kD9k+At82fn8skyLY46hAzM7scuBzgkd3DWWksQArQC0RNbxpHOYVYVPb+aITdcpeeqcTiVmwFUTbJ0QRNlOMKscQHxNrIpfhI8QGxzlkV0d5nkl0bxLzj7IoAM7NXME2RfsuJfu88VO0AwLn7vtybCzBRqaCq/C6UAYXIWTDSJjn7qnREEzWRNgaRNnLRxEekvzEVkc6vKER6f4lGpL/ZYG+HD5AZMAFm9mImc45n+YNH+Fbg/G0ve0hDzAxY7drq5IKm9yda/5PiDSZSmVwkQZPCYyLSJi5Fgo5I54GK3HxriPR3EI1o74mJlsyANcbMLgJ+HPhWdz+47UvXAr9mZj/LZMJxAfDHx4tXzDUCTJX9CSRqIvU/RRJg0cRHpItipM1RNJGQG/pY51ckIr3HRCPa+0ySwJYL4ul/7ipt6N8KPAM4y8xuAV7F5Hq4D3i3Te+6H3D3H3T3j5rZNcDHmEoTf+h4DogwZcD2lbYCTDnMVSVqogiQFB86VJvDSBfwFAkpEiDW3200Ir0fJPmemMR9P9zpOWC7gdIF8flHePqok6Hd/bXAa0/kZ0xzwBpnwAIJGoDONIla1R9tlA2icqMR6aIY5felJOoFrCW58dYR6f0giUW+fydRXRCzBHGPU8zZ3w+SuM1jyrI0cXrAVETbIC/7ZjY3nDoi/d0m8bAlf+9KkkSP49QleK+JLcBwVgUCTEFmVHTkplNHbrgSFfl3myRJ8iCKm8nB7k9PePaA7XnMnEU5bqvYicXMDSeQm6MkHtEyoUmSJEmSfCnZA7bHMYN9jTNg7rmLS5IkSZIkSZKdZhrEnAJsT/OJ+++84zkffPv9wB27vZZkT3IWeW4kRybPjeRY5PmRHI08N5KjsRfPja/c7QWcDCnA9jju/mVmdr27X7jba0n2HnluJEcjz43kWOT5kRyNPDeSo5HnRit8KUoQVe7oSZIkSZIkSZIkyWGEzoAlSZIkSZIkSXJ6kD1gcTiw2wtI9ix5biRHI8+N5Fjk+ZEcjTw3kqOR50YLDKqd/qOYzZfAaz9JkiRJkiRJkr3Nw7qz/GtWv+OEv+9PDr75hkg9eKdDBixJkiRJkiRJkuBMFhynfwYsBViSJEmSJEmSJHuCZegBC+uCaGYXmdlfmNmNZnbFbq8n2XnM7E1mdruZfWTbc48xs3eb2Sfm/z96ft7M7D/N58v/MLNv3L2VJ2rM7Hwze5+ZfczMPmpmPzo/n+fHkmNmq2b2x2b2Z/O58W/n559oZh+cz4G3mdnK/Py++fGN89efsKv/gESOmXVm9qdm9l/mx3luJACY2afN7MNm9iEzu35+Lq8rjalWT/gjGiEFmJl1wC8A3wY8CXi+mT1pd1eV7AK/Alx02HNXAO9x9wuA98yPYTpXLpg/Lgdev0NrTHaHAfjX7v4k4JuBH5rfI/L8SNaBZ7r7NwBPAS4ys28Gfhr4OXf/G8CdwGXz6y8D7pyf/7n5dcnpzY8CH9/2OM+NZDv/wN2fsq3fKK8rDZkKEE/8v2iEFGDA04Ab3f0md98ArgYu3uU1JTuMu/8B8MXDnr4YuGr+/Crgu7Y9/2af+ABwppmdsyMLTXYcd7/N3f9k/vxeps3UueT5sfTMv+P75oeL+cOBZwJvn58//NzYOmfeDjzLzGxnVpvsNGZ2HvAdwBvmx0aeG8mxyetKY1KA7V3OBW7e9viW+bkkOdvdb5s//xxw9vx5njNLylwW9FTgg+T5kfBAidmHgNuBdwOfBO5y92F+yfbf/wPnxvz1u4HH7uiCk53kPwI/Dg/s6B5LnhvJgzjwe2Z2g5ldPj+X15WmTDmwE/2IRlQBliTHxacZC6d/J2dyVMzs4cBvAP/C3e/Z/rU8P5YXdx/d/SnAeUwVFV+zuytK9gJm9p3A7e5+w26vJdmz/K/u/o1M5YU/ZGbfsv2LeV05dZyd7wE7Wh/fYa/5B3Pv39bHmpl91/y1XzGzT2372lOO9zOjCrBbgfO3PT5vfi5JPr+V4p//f/v8fJ4zS4aZLZjE11vc/Tfnp/P8SB7A3e8C3gf8XabyoC1n4O2//wfOjfnrjwK+sLMrTXaIpwPPNbNPM7U2PBN4HXluJDPufuv8/9uB32K6gZPXlabsSg/Y0fr4HlyV+/vm3r+nML03HAR+b9tLfmzr6+7+oeP9wKgC7DrggtmZaAW4FLh2l9eU7A2uBV40f/4i4B3bnv/+2ZXom4G7t5UMJKcZcx/GG4GPu/vPbvtSnh9Ljpl9mZmdOX++H3gOU4/g+4DnzS87/NzYOmeeB7x3vsudnGa4+8vd/Tx3fwLTvuK97v595LmRAGb2MDN7xNbnwD8EPkJeV5rjjCf8cYocrY/vaDwP+G/ufvBkf2DIOWDuPpjZS4F3AR3wJnf/6C4vK9lhzOytwDOAs8zsFuBVwJXANWZ2GfAZ4JL55e8Evh24kemuxQ/s+IKTneTpwAuBD8+9PgA/QZ4fCZwDXDW76RbgGnf/L2b2MeBqM/sp4E+ZBDzz/3/VzG5kMv25dDcWnewqLyPPjWTq7fqt2WelB37N3X/XzK4jryvN2HJB3GGO1sd3NC4Ffvaw515rZq9kzqC5+/qxAljerEmSJEmSJEmSZLfZ1z3Kzznj753w933mvt/9DHDHtqcOuPuBrQdm9vvA447wra8ArnL3M7e99k53/5I+sPlr5wD/A3i8u29ue+5zwApwAPiku7/mWOsNmQFLkiRJkiRJkuR0w0+2pPCObbPZvjSq+7OP9jUz+7yZnePutx3Wx3ckLgF+a0t8zbG3smfrZvbLwL853mKj9oAlSZIkSZIkSXIa4ezKHLCj9fEdiecDb93+xDYTFmPqH/vI8X5gCrAkSZIkSZIkSZaVK4HnmNkngGfPjzGzC83sDVsvmueKng/8v4d9/1vM7MPAh4GzgJ863g/MHrAkSZIkSZIkSXadle6R/uVnHLWS8Kjcet/7bjhWCeJeI3vAkiRJkiRJkiTZAzj11G3l9zxZgpgkSZIcFTM708z++fz5483s7bu9piRJkuT0xJkk2Il+RCMFWJIkSXIszgT+OYC7/093f96xX54kSZIkJ4tTfTzhj2hkCWKSJElyLK4EvnoeaP0J4Gvd/evM7MVMbk8PAy4A/gPTDJQXAuvAt7v7F83sq4FfAL6MaRjpP3X3P9/pf0SSJEkSg4gZrRMlM2BJkiTJsbiCaajkU4AfO+xrXwf878DfAV4LHHT3pwJ/BHz//JoDwA+7+99mmo3y/+zEopMkSZKITHPATvQjGpkBS5IkSU6W97n7vcC9ZnY38Dvz8x8GnmxmDwf+HvDr03gUAPbt/DKTJEmSCDhQ/fTPgKUAS5IkSU6W9W2f122PK9P1pQB3zdmzJEmSJDkOniWISZIkydJzL/CIk/lGd78H+JSZfTeATXxDy8UlSZIkpxEO7uMJf0QjM2BJkiTJUXH3L5jZ/2dmHwE+fhIhvg94vZn9JLAArgb+rOUakyRJktODKf91+mfAzN13ew1JkiRJkiRJkiw5XdnvD9v3xBP+vnvXPn6Du18oWJKEzIAlSZIkSZIkSbIH8JCuhidKCrAkSZIkSZIkSfYEni6ISZIkSZIkSZIkO8FyuCCmAEuSJEmSJEmSZNdxCOlqeKKkAEuSJEmSJEmSZA/gWYKYJEmSJEmSJEmyUyxDCWIOYk6SJEmSJEmSJNkhMgOWJEmSJEmSJMnu4+mCmCRJkiRJkiRJskOkC2KSJEmSJEmSJMmOkC6ISZIkSZIkSZIkO4ZDZsCSJEmSJEmSJEl2huwBS5IkSZIkSZIk2RGyByxJkiRJkiRJkmQHSQGWJEmSJEmSJEmyM2QJYpIkSZIkSZIkyU6QJYhJkiRJkiRJkiQ7SAqwJEmSJEmSJEmSncF9t1cgJwVYkiRJkiRJkiR7AMdJAZYkSZIkSZIkSbITvAuGs07i++5ovhIh5kuQ5kuSJEmSJEmSJNkLlN1eQJIkSZIkSZIkybKQAixJkiRJkiRJkmSHSAGWJEmSJEmSJEmyQ6QAS5IkSZIkSZIk2SFSgCVJkiRJkiRJkuwQ/z8oaw6lxKrYyAAAAABJRU5ErkJggg==" }, "metadata": { "needs_background": "light" } } ], "metadata": {} }, { "cell_type": "markdown", "source": [ "This concludes a first simulation in phiflow. It's not overly complex, but because of that it's a good starting point for evaluating and comparing different physics-based deep learning approaches in the next chapter. But before that, we'll target a more complex simulation type in the next section." ], "metadata": {} }, { "cell_type": "markdown", "source": [ "## Next steps\n", "\n", "Some things to try based on this simulation setup:\n", "\n", "- Feel free to experiment - the setup above is very simple, you can change the simulation parameters, or the initialization. E.g., you can use a noise field via `Noise()` to get more chaotic results (cf. the comment in the `velocity` cell above).\n", "\n", "- A bit more complicated: extend the simulation to 2D (or higher). This will require changes throughout, but all operators above support higher dimensions. Before trying this, you probably will want to check out the next example, which covers a 2D Navier-Stokes case." ], "metadata": {} } ], "metadata": { "kernelspec": { "display_name": "Python 3", "language": "python", "name": "python3" }, "language_info": { "codemirror_mode": { "name": "ipython", "version": 3 }, "file_extension": ".py", "mimetype": "text/x-python", "name": "python", "nbconvert_exporter": "python", "pygments_lexer": "ipython3", "version": "3.8.5" } }, "nbformat": 4, "nbformat_minor": 4 }