{ "cells": [ { "attachments": {}, "cell_type": "markdown", "metadata": {}, "source": [ "# Time-varying depths\n" ] }, { "attachments": {}, "cell_type": "markdown", "metadata": {}, "source": [ "Some hydrodynamic models (such as SWASH) have time-evolving depth dimensions, for example because they follow the waves on the free surface. Parcels can work with these types of models, but it is a bit involved to set up. That is why we explain here how to run Parcels on `FieldSets` with time-evolving depth dimensions\n" ] }, { "cell_type": "code", "execution_count": 1, "metadata": {}, "outputs": [], "source": [ "from datetime import timedelta as delta\n", "from os import path\n", "\n", "import matplotlib.pyplot as plt\n", "import numpy as np\n", "import xarray as xr\n", "\n", "from parcels import (\n", " AdvectionRK4,\n", " FieldSet,\n", " JITParticle,\n", " ParticleFile,\n", " ParticleSet,\n", " download_example_dataset,\n", ")" ] }, { "attachments": {}, "cell_type": "markdown", "metadata": {}, "source": [ "Here, we use sample data from the SWASH model. We first set the `filenames` and `variables`\n" ] }, { "cell_type": "code", "execution_count": 2, "metadata": {}, "outputs": [], "source": [ "example_dataset_folder = download_example_dataset(\"SWASH_data\")\n", "filenames = f\"{example_dataset_folder}/field_*.nc\"\n", "variables = {\n", " \"U\": \"cross-shore velocity\",\n", " \"V\": \"along-shore velocity\",\n", " \"depth_u\": \"time varying depth_u\",\n", "}" ] }, { "attachments": {}, "cell_type": "markdown", "metadata": {}, "source": [ "Now, the first key step when reading time-evolving depth dimensions is that we specify `depth` as `'not_yet_set'` in the `dimensions` dictionary\n" ] }, { "cell_type": "code", "execution_count": 3, "metadata": {}, "outputs": [], "source": [ "dimensions = {\n", " \"U\": {\"lon\": \"x\", \"lat\": \"y\", \"depth\": \"not_yet_set\", \"time\": \"t\"},\n", " \"V\": {\"lon\": \"x\", \"lat\": \"y\", \"depth\": \"not_yet_set\", \"time\": \"t\"},\n", " \"depth_u\": {\"lon\": \"x\", \"lat\": \"y\", \"depth\": \"not_yet_set\", \"time\": \"t\"},\n", "}" ] }, { "attachments": {}, "cell_type": "markdown", "metadata": {}, "source": [ "Then, _after_ we create the `FieldSet` object, we set the `depth` dimension of the relevant `Fields` to `fieldset.depth_u` and `fieldset.depth_w`, using the `Field.set_depth_from_field()` method\n" ] }, { "cell_type": "code", "execution_count": 4, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "WARNING: Flipping lat data from North-South to South-North. Note that this may lead to wrong sign for meridional velocity, so tread very carefully\n" ] } ], "source": [ "fieldset = FieldSet.from_netcdf(\n", " filenames, variables, dimensions, mesh=\"flat\", allow_time_extrapolation=True\n", ")\n", "fieldset.U.set_depth_from_field(fieldset.depth_u)\n", "fieldset.V.set_depth_from_field(fieldset.depth_u)" ] }, { "attachments": {}, "cell_type": "markdown", "metadata": {}, "source": [ "Now, we can create a ParticleSet, run those and plot them\n" ] }, { "cell_type": "code", "execution_count": 5, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "INFO: Output files are stored in SwashParticles.zarr.\n", "100%|██████████| 0.25/0.25 [00:00<00:00, 1.61it/s] \n" ] } ], "source": [ "pset = ParticleSet(fieldset, JITParticle, lon=9.5, lat=12.5, depth=-0.1)\n", "\n", "pfile = pset.ParticleFile(\"SwashParticles\", outputdt=delta(seconds=0.05))\n", "\n", "pset.execute(AdvectionRK4, dt=delta(seconds=0.005), output_file=pfile)" ] }, { "cell_type": "code", "execution_count": 6, "metadata": {}, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAkcAAAGsCAYAAADXKAMrAAAAOXRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjguMCwgaHR0cHM6Ly9tYXRwbG90bGliLm9yZy81sbWrAAAACXBIWXMAAA9hAAAPYQGoP6dpAABS8ElEQVR4nO3deVhU5/k+8Hv2AYRBFllUBlAD4g64QGNME4OaJtGYBExSql1STZsYzCYurek3bcRsan9xSazZ2ia4EMzWGEmjxAWIGMAFFwTcQQRh2HSGmXl/f1gnTkBklOGw3J/rmj94z3Pe85w5zjW3Z86ckQkhBIiIiIgIACCXugEiIiKizoThiIiIiOgaDEdERERE12A4IiIiIroGwxERERHRNRiOiIiIiK7BcERERER0DYYjIiIiomswHBERERFdg+GIiIiI6BoMRx3knXfewZ133gkPDw/IZDLU1NTccJ2lS5di9OjRcHd3R58+fTBt2jQcPXr0uvWzZ8+GTCbDihUr7MaLi4vx4IMPwtfXFx4eHoiPj8f58+ftan744Qfcc8898PT0hLe3N37/+9+jvr7eoX387rvvcP/99yMwMBAymQxbtmxxaH0iIqLOgOGoHd155514//33W1zW2NiIyZMnY+HChW2eLzMzE3/84x+RnZ2NjIwMmM1mxMXFoaGhoVntli1bkJOTg8DAQLvxhoYGxMXFQSaT4dtvv8Xu3bthMplw//33w2q1AgDOnTuHiRMnYuDAgcjJycHWrVtx6NAhzJo1q829Xt3WiBEj8NZbbzm0HhERUaciqN1MmDBBvPfee63WbN++XQAQ1dXVDs9fUVEhAIjMzEy78TNnzoi+ffuKgwcPCr1eL5YvX25b9vXXXwu5XC4MBoNt7OLFiwKAyMjIEEII8fbbb4s+ffoIi8Viq8nLyxMARFFRkW3s0KFDYsqUKcLNzU306dNH/PKXvxQXLlxosVcAIj093eF9JCIikhrPHHUhBoMBAODl5WUbs1qtSExMxAsvvIAhQ4Y0W8doNEImk0Gj0djGtFot5HI5du3aZatRq9WQy3/85+Di4gIAtpqysjJMmDABI0eORG5uLrZu3Yrz588jPj6+/XeUiIhIQgxHXYQQAs8++yxuv/12DB061Da+bNkyKJVKzJ07t8X1xo0bBzc3N8yfPx+NjY1oaGjACy+8AKvVirKyMgDAXXfdhfLycrz22mswmUyorq62ffx3tWbNmjWIjIzEK6+8gvDwcIwaNQrvvvsutm/fjmPHjjl574mIiDoOw9EteOWVV9CrVy/bY+fOnZgzZ06zsfbw1FNPYf/+/fj4449tY/v27cPKlSvx/vvvQyaTtbier68vNm3ahM8//xy9evWCTqeDwWBAZGQkFAoFAGDIkCH44IMP8MYbb8DV1RX+/v4IDQ2Fn5+frWbfvn3Yvn273b6Fh4cDuHLBNxERUXehlLqBrmzOnDl2Hys9/vjjeOihhzB9+nTbWN++fW95O08//TQ+++wzfPfdd+jXr59tfOfOnaioqEBQUJBtzGKx4LnnnsOKFStw4sQJAEBcXByKi4tRWVkJpVIJT09P+Pv7IyQkxLbeY489hsceewznz5+Hm5sbZDIZ3nzzTVuN1WrF/fffj2XLljXrLyAg4Jb3kYiIqLNgOLoFXl5edtf/uLi4oE+fPhg4cGC7zC+EwNNPP4309HTs2LHDLswAQGJiIiZOnGg3NmnSJCQmJuLXv/51s/l8fHwAAN9++y0qKirwwAMPNKvx8/MDALz77rvQarW45557AACRkZFIS0tDcHAwlEr+syEiou6L73IdpLy8HOXl5Th+/DgA4MCBA3B3d0dQUJAtYN1999148MEH8dRTTwEA/vjHP+Kjjz7Cp59+Cnd3d5SXlwMAdDodXFxc4O3tDW9vb7vtqFQq+Pv7IywszDb23nvvYfDgwfD19UVWVhaeeeYZzJs3z67mrbfeQmxsLHr16oWMjAy88MILSElJgaenp62XdevW4dFHH8ULL7wAHx8fHD9+HKmpqVi3bh0UCgXq6+tt+wcApaWlyM/Ph5eXl93ZLSIiok5N6q/LdSetfZV/yZIlAkCzx7X1er1eLFmyxPZ3S/U/XeenfvpVfiGEmD9/vvDz8xMqlUoMGjRIvPHGG8JqtdrVJCYmCi8vL6FWq8Xw4cPFhx9+2GzuY8eOiQcffFB4enoKFxcXER4eLpKSkmxzXb1NwU8fM2fObO1pIyIi6lRkQgjR8ZGMiIiIqHPit9WIiIiIrsFwRERERHQNXpDtIKvVinPnzsHd3f269xYiIiKizkUIgbq6OgQGBtr9IkRLGI4cdO7cOfTv31/qNoiIiOgmnD592u6egS1hOHKQu7s7gCtProeHh8TdEBERUVvU1taif//+tvfx1jAcOejqR2keHh4MR0RERF1MWy6J4QXZRERERNdgOCIiIiK6BsMRERER0TUYjoiIiIiuwXBEREREdA2GIyIiIqJrMBwRERERXYPhiIiIiOgaDEdERERE17ipcLR69WqEhIRAq9UiKioKO3fubLU+MzMTUVFR0Gq1CA0Nxdq1a5vVpKWlISIiAhqNBhEREUhPT7db/t133+H+++9HYGAgZDIZtmzZ0mwOIQReeuklBAYGwsXFBXfeeScOHTpkV2M0GvH000/Dx8cHbm5ueOCBB3DmzBnHnwQiIiLqlhwORxs2bEBSUhIWLVqEvLw8jB8/HlOmTMGpU6darC8tLcW9996L8ePHIy8vDwsXLsTcuXORlpZmq8nKykJCQgISExNRUFCAxMRExMfHIycnx1bT0NCAESNG4K233rpub6+++irefPNNvPXWW9i7dy/8/f1xzz33oK6uzlaTlJSE9PR0pKamYteuXaivr8d9990Hi8Xi6FPR7soMl7CnuBJlhktSt0JERNRzCQeNGTNGzJkzx24sPDxcJCcnt1j/4osvivDwcLux2bNni3Hjxtn+jo+PF5MnT7armTRpkpgxY0aLcwIQ6enpdmNWq1X4+/uLlJQU29jly5eFTqcTa9euFUIIUVNTI1QqlUhNTbXVnD17VsjlcrF169br7LE9g8EgAAiDwdCm+rZK/f6kCEn+QujnfyFCkr8Qqd+fbNf5iYiIejJH3r8dOnNkMpmwb98+xMXF2Y3HxcVhz549La6TlZXVrH7SpEnIzc1FU1NTqzXXm7MlpaWlKC8vt5tHo9FgwoQJtnn27duHpqYmu5rAwEAMHTr0utsyGo2ora21e7S3MsMlLPjkAKziyt9WASz45ADPIBEREUnAoXBUWVkJi8UCPz8/u3E/Pz+Ul5e3uE55eXmL9WazGZWVla3WXG/O623n6nrXm6e8vBxqtRq9e/du87aWLl0KnU5ne/Tv37/NPbVVaWWDLRhdZRXA0x/l4b3dpThSXgvrTwuIiIjIKZQ3s5JMJrP7WwjRbOxG9T8dd3TO9urtRjULFizAs88+a/u7tra23QNSiI8b5DI0C0i5J6uRe7IaAODlpsa4UC/EhHojZoA3Bvj2uqnnh4iIiFrnUDjy8fGBQqFodpaloqKi2Rmbq/z9/VusVyqV8Pb2brXmenNebzvAlbNDAQEBLc7j7+8Pk8mE6upqu7NHFRUViI2NbXFejUYDjUbT5j5uRoDOBUunD8PCTw7CIgTkMuDJCQPgplUiq7gKuSeqcbHBhP8cKMd/Dlx5nnzdNRgX6m0LS8HergxLRERE7cChcKRWqxEVFYWMjAw8+OCDtvGMjAxMnTq1xXViYmLw+eef241t27YN0dHRUKlUtpqMjAzMmzfPruZ6gaUlISEh8Pf3R0ZGBkaNGgXgyjVSmZmZWLZsGQAgKioKKpUKGRkZiI+PBwCUlZXh4MGDePXVV9u8LWdIGB2EO27zxYnKRgT7uCJA5wIA+MOdA2EyW7H/TA2yiquQVVKFfSercaHOiM8LzuHzgnMAAH8PLWIG/BiW+nu5Srk7REREXZbDH6s9++yzSExMRHR0NGJiYvDOO+/g1KlTmDNnDoArH0OdPXsWH374IQBgzpw5eOutt/Dss8/iiSeeQFZWFtavX4+PP/7YNuczzzyDO+64A8uWLcPUqVPx6aef4ptvvsGuXbtsNfX19Th+/Ljt79LSUuTn58PLywtBQUGQyWRISkrCK6+8gkGDBmHQoEF45ZVX4OrqisceewwAoNPp8Nvf/hbPPfccvL294eXlheeffx7Dhg3DxIkTb+4ZbEcBOhdbKLqWWilHdLAXooO98PTdg3C5yYL80z+GpfxTNSivvYz0vLNIzzsLAOjr6WIXlgI9m89LRERELbiZr8OtWrVK6PV6oVarRWRkpMjMzLQtmzlzppgwYYJd/Y4dO8SoUaOEWq0WwcHBYs2aNc3m3LRpkwgLCxMqlUqEh4eLtLQ0u+Xbt28XAJo9Zs6caauxWq1iyZIlwt/fX2g0GnHHHXeIAwcO2M1z6dIl8dRTTwkvLy/h4uIi7rvvPnHq1Kk277uzvsp/KxqNZrGr6IJ4besRMX31bjFgwZdCP/8Lu8cdr34r5m8uEFvyzojzhktSt0xERNShHHn/lgkh+DUoB9TW1kKn08FgMMDDw0PqdlrUYDQj92S17czSgTM1zS72DvV1s51VGhfqDZ9ezr2uioiISEqOvH8zHDmoK4Sjn6q73IS9Jy7awtKhc7X46VG/za+XLSyNDfFGbze1NM0SERE5AcORE3XFcPRThsYm5JRWYU9xFbJLqnCkvM5uuUwGhPt72MLSmBAv6FxUEnVLRER06xiOnKg7hKOfqqo3Iqf0xzNLxyvq7ZbLZcCQQJ3tAu/RIV7opbmpW2QRERFJguHIibpjOPqpirrLyC65EpayS6pQWtlgt1whl2FY3x/DUnRwb7iqGZaIiKjzYjhyop4Qjn6q3HAZWSWVtjNLpy/a/+abSiHDiH6etrAUqe8NrUohUbdERETNMRw5UU8MRz91prrRFpSyi6twznDZbrlaKUdkkCdiQn0QM8AbI/rroFEyLBERkXQYjpyI4cieEAKnLv4YlrKKq1BRZ7Sr0arkiNZ72W4bMLyfDiqFQ795TEREdEsYjpyI4ah1QgiUVDbYwlJOSRUq6012Na5qBUYHe9k+hhsS6AElwxIRETkRw5ETMRw5RgiBoor6K2GpuArZpVWoaWyyq3HXKDEm5MczSxEBHpDL+SO6RETUfhiOnIjh6NZYrQJHyutsH8HllFah7rLZrkbnosLY/4WlmAHeuK2PO8MSERHdEoYjJ2I4al8Wq0DhuVrbt+H2nqhGvdE+LHm5qTEu1Mt2U8oBvr0gkzEsERFR2zEcORHDkXOZLVYcOGuwnVnKPVGNS00Wuxpfdw3GhXrbwlKwtyvDEhERtYrhyIkYjjqWyWzF/jM1tgu8952shtFstavx99DaLu6OGeCN/l6uEnVLRESdFcOREzEcSetykwX5p38MS/mnamCy2Ielvp4udmEp0NNFom6JiKizYDhyIoajzuWSyYIfTlXbwlLB6RqYrfb/pPXerragFBPqjT4eWom6JSIiqTAcORHDUefWYDQj9+SPYenAmRr8JCsh1NfNFpbGhXrDp5dGmmaJiKjDMBw5EcNR11J3uQl7T1y0haVD52rx03/xt/n1soWlsSHe6O2mlqZZIiJyGoYjJ2I46toMjU3IKf3xp06OlNfZLZfJgHB/D1tYGhPiBZ2LSqJuiYiovTAcORHDUfdyscGEnJIfw1JRRb3dcrkMGBKos12vNDrEC700Som6JSKim8Vw5EQMR93bhTojsv8XlrKLq1BS2WC3XCGXYVjfH8NSdHBvuKoZloiIOjuGIydiOOpZyg2Xr4Sl/12zdOpio91ylUKGEf08bT91EhnUG1qVQqJuiYjoehiOnIjhqGc7U92I7JKL2FNcieziKpwzXLZbrlbKERnkiZhQH8QM8MaI/jpolAxLRERSYzhyIoYjukoIgVMXG21nlbKKq1BRZ7Sr0arkiNZ72W4bMLyfDiqFXKKOiYh6LoYjJ2I4ousRQqCkssEWlnJKqlBZb7KrcVUrMDrYy3bN0pBADygZloiInI7hyIkYjqithBAoqqi/EpaKq5BdWoWaxia7GneNEmNCfjyzFBHgAbmcP6JLRNTeGI6ciOGIbpbVKnCkvM72EVxOaRXqLpvtanQuKoz9X1iKGeCN2/q4MywREbUDhiMnYjii9mKxChSeq0VWSSWyiquw90Q16o32YcnLTY1xoV62m1IO8O0FmYxhiYjIUQxHTsRwRM5itlhx4KzBdmYp90Q1LjVZ7Gp83TUYF+ptC0vB3q4MS0REbcBw5EQMR9RRTGYr9p+psV3gve9kNYxmq12Nv4fWdnF3zABv9PdylahbIqLOjeHIiRiOSCqXmyzIP/1jWMo/VQOTxT4s9fV0sQtLgZ4uEnVLRNS5MBw5EcMRdRaXTBb8cKraFpYKTtfAbLV/Oeu9XW1BKSbUG308tBJ1S0QkLYYjJ2I4os6qwWhG7skfw9KBMzX4SVZCqK+bLSyNC/WGTy+NNM0SEXUwhiMnYjiirqLuchP2nrhoC0uHztXip6/22/x62cLS2BBv9HZTS9MsEZGTMRw5EcMRdVWGxibklP74UydHyuvslstkQLi/hy0sjQnxgs5FJVG3RETti+HIiRiOqLu42GBCTsmPYamoot5uuVwGDAnU2a5XGh3ihV4apUTdEhHdGoYjJ2I4ou7qQp0R2f8LS9nFVSipbLBbrpDLMKzvj2EpOrg3XNUMS0TUNTAcORHDEfUU5YbLV8LS/65ZOnWx0W65SiHDiH6etp86iQzqDa1KIVG3REStYzhyIoYj6qnOVDciu+Qi9hRXIru4CucMl+2Wq5VyRAZ5IibUBzEDvDGivw4aJcMSEXUODEdOxHBEBAghcOpio+2sUlZxFSrqjHY1WpUc0Xov220DhvfTQaWQS9QxEfV0DEdOxHBE1JwQAiWVDbawlFNShcp6k12Nm1qB6GAv2zVLQwI9oGRYIqIOwnDkRAxHRDcmhEBRRf2VsFRchezSKtQ0NtnVuGuUGBPy45mliAAPyOX8EV0icg6GIydiOCJynNUqcKS8zvYRXE5pFeoum+1qdC4qjP1fWIoZ4I3b+rgzLBFRu2E4ciKGI6JbZ7EKFJ6rRVZJJbKKq7D3RDXqjfZhyctNjXGhXrabUg7w7QWZjGGJiG4Ow5ETMRwRtT+zxYoDZw22M0u5J6pxqcliV+PrrsG4UG9bWAr2drWFpTLDJZRWNiDExw0BOhcpdoGIOjmGIydiOCJyPpPZiv1namwXeO87WQ2j2WpX4++hRcwAbyhkMnySdwZWceWu3kunD0PC6CCJOieizorhyIkYjog6ntFsQd6pH8NS/qkamCzWFmsVMhl2Jf+cZ5CIyI4j79+89z8RdXoapQLjQq98q20egMtNFuw7WY1NuaexJf+cXa1FCJyobGQ4IqKbxpuMEFGXo1Up8LOBPpg/JRwtfaHt+9Kqjm+KiLoNhiMi6rICdC5YOn0YFP+7MPtqTlr+TRHe3HYUvGqAiG4GP1Yjoi4tYXQQ7rjNFycqG6H3dsHmfWfxZsYx/P3b47hQb8Jfpw2FgvdLIiIHMBwRUZcXoHOxXWM09+5B8O6lxuItB/Hx96dQ3WDCihkjoVXxR3CJqG34sRoRdTuPj9Vj9WORUCvk2HqoHLPe+x51l5tuvCIRERiOiKibmjIsAO//ZjR6aZTILrmIGe9k40KdUeq2iKgLYDgiom4rdoAPUn8/Dj691Dh0rhYPr92DU1WNUrdFRJ0cwxERdWtD++qweU4s+nu54GRVIx5auweF52qlbouIOjGGIyLq9oJ93JA2Jxbh/u64UGdEwttZyCnhvZCIqGUMR0TUI/Tx0GLD7BiMCfZCndGMxHe/x7ZD5VK3RUSdEMMREfUYOhcVPvztGNwT4QeT2Yo5/9qHDXtPSd0WEXUyDEdE1KNoVQqseTwS8dH9YBXA/LQDWL3jOO+mTUQ2NxWOVq9ejZCQEGi1WkRFRWHnzp2t1mdmZiIqKgparRahoaFYu3Zts5q0tDRERERAo9EgIiIC6enpDm/3/PnzmDVrFgIDA+Hq6orJkyejqKjIrqa8vByJiYnw9/eHm5sbIiMjsXnz5pt4Foioq1Iq5Fj20HA8eecAAMCrW4/i5S8Ow2plQCKimwhHGzZsQFJSEhYtWoS8vDyMHz8eU6ZMwalTLZ+aLi0txb333ovx48cjLy8PCxcuxNy5c5GWlmarycrKQkJCAhITE1FQUIDExETEx8cjJyenzdsVQmDatGkoKSnBp59+iry8POj1ekycOBENDQ22eRITE3H06FF89tlnOHDgAKZPn46EhATk5eU5+lQQURcmk8kwf3I4Fv9iMADg3d2leHZjPposVok7IyLJCQeNGTNGzJkzx24sPDxcJCcnt1j/4osvivDwcLux2bNni3Hjxtn+jo+PF5MnT7armTRpkpgxY0abt3v06FEBQBw8eNC23Gw2Cy8vL7Fu3TrbmJubm/jwww/t5vHy8hL/+Mc/rrvP1zIYDAKAMBgMbaonos7vkx9OiwELvhT6+V+IX63PEQ3GJqlbIqJ25sj7t0NnjkwmE/bt24e4uDi78bi4OOzZs6fFdbKysprVT5o0Cbm5uWhqamq15uqcbdmu0Xjlzrdarda2XKFQQK1WY9euXbax22+/HRs2bMDFixdhtVqRmpoKo9GIO++8s8X+jUYjamtr7R5E1L08OKof1s2MhlYlR+axC3hsXQ6qG0xSt0VEEnEoHFVWVsJiscDPz89u3M/PD+XlLX8ltry8vMV6s9mMysrKVmuuztmW7YaHh0Ov12PBggWorq6GyWRCSkoKysvLUVZWZltnw4YNMJvN8Pb2hkajwezZs5Geno4BAwa02P/SpUuh0+lsj/79+9/oaSKiLujnYX3w79+Ng85FhfzTNXjk7Sycq7kkdVtEJIGbuiBbJpPZ/S2EaDZ2o/qfjrdlztZqVCoV0tLScOzYMXh5ecHV1RU7duzAlClToFD8+GvcixcvRnV1Nb755hvk5ubi2WefxSOPPIIDBw602PuCBQtgMBhsj9OnT193P4moa4vS98bmOTEI0GlxvKIeD63Zg+MVdVK3RUQdTOlIsY+PDxQKRbOzRBUVFc3O6lzl7+/fYr1SqYS3t3erNVfnbOt2o6KikJ+fD4PBAJPJBF9fX4wdOxbR0dEAgOLiYrz11ls4ePAghgwZAgAYMWIEdu7ciVWrVrX4LTqNRgONRnPD54aIuodBfu7Y/GQsfrU+B8UXGvDw2iy8N2s0RgX1lro1IuogDp05UqvViIqKQkZGht14RkYGYmNjW1wnJiamWf22bdsQHR0NlUrVas3VOR3drk6ng6+vL4qKipCbm4upU6cCABobr/zgpFxuv9sKhQJWK7+hQkRX9PV0waY5sRjR3xM1jU14bF0OMo9dkLotIuoojl7tnZqaKlQqlVi/fr0oLCwUSUlJws3NTZw4cUIIIURycrJITEy01ZeUlAhXV1cxb948UVhYKNavXy9UKpXYvHmzrWb37t1CoVCIlJQUcfjwYZGSkiKUSqXIzs5u83aFEGLjxo1i+/btori4WGzZskXo9Xoxffp023KTySQGDhwoxo8fL3JycsTx48fF66+/LmQymfjyyy/btP/8thpRz1F/uUkkrs8R+vlfiAELvhRb8s5I3RIR3SRH3r8dDkdCCLFq1Sqh1+uFWq0WkZGRIjMz07Zs5syZYsKECXb1O3bsEKNGjRJqtVoEBweLNWvWNJtz06ZNIiwsTKhUKhEeHi7S0tIc2q4QQqxcuVL069dPqFQqERQUJBYvXiyMRqNdzbFjx8T06dNFnz59hKurqxg+fHizr/a3huGIqGcxNlnE0x/9IPTzvxD6+V+Id3eVSN0SEd0ER96/ZULwnvmOqK2thU6ng8FggIeHh9TtEFEHsFoF/u+LQry/5wQA4I8/H4Dn48Ja/SIKEXUujrx/87fViIhuQC6XYcn9EXhhUhgAYNX2Yiz45ADMvJs2UbfEcERE1AYymQx//PlALJ0+DHIZkLr3NP7w7x9wuckidWtE1M4YjoiIHPDomCCsfjwKaqUc2wrPY+a736P2cpPUbRFRO2I4IiJy0OSh/vjg12PgrlEip/QiEt7ORkXdZanbIqJ2wnBERHQTYgZ4I3X2OPj00uBwWS0eXpOFk1UNUrdFRO2A4YiI6CYNCdQh7ckYBHm54tTFRjy0JgsHzxqkbouIbhHDERHRLdB7u2HzkzEYHOCBynojHn0nG1nFVVK3RUS3gOGIiOgW9XHXYsPscRgb4oU6oxkz3/seWw+WSd0WEd0khiMionbgoVXhg9+MQVyEH0xmK/7w7x/w8fenpG6LiG4CwxERUTvRqhRY/XgkZozuD6sAFnxyAG99WwT+EAFR18JwRETUjpQKOZZOH4anfj4QAPD6tmP4y+eFsFoZkIi6CoYjIqJ2JpPJ8PykMCy5PwIA8P6eE0jakA+TmT83QtQVMBwRETnJr38WgpUzRkIpl+GzgnP47Qd70WA0S90WEd0AwxERkRNNHdkX/5gZDReVAjuLKvHYP3JwscEkdVtE1AqGIyIiJ7szrA8+emIsPF1VKDhdg4fX7sHZmktSt0VE18FwRETUAUYF9cbmOTEI1GlRcqEBD6/Zg6LzdVK3RUQtYDgiIuogA/u4Y/OTsRjYpxfKDJfxyNtZ2HeyWuq2iOgnGI6IiDpQoKcLNs2OwaggT9Q0NuGX/8jB9qMVUrdFRNdgOCIi6mC93dT49+/GYsJtvrjUZMETH+RiS95Zqdsiov9hOCIikoCrWol/zIzGtJGBMFsFkjbkY/2uUqnbIiIwHBERSUalkOPN+JH4zc9CAAAvf1GIZVuP8OdGiCTGcEREJCG5XIY/3TcYL04OAwCs2VGM+Wn7YbbwbtpEUmE4IiKSmEwmwx/uHIhlDw2DXAZszD2DJ//9Ay43WaRujahHYjgiIuokEkYHYc0vo6BWypFReB6/Wv89DJeapG6LqMdhOCIi6kQmDfHHP38zBu4aJb4/cREJb2ehovay1G0R9SgMR0REnczYUG9smB0DX3cNjpTX4aG1e3CiskHqtoh6DIYjIqJOKCLQA2lzYqH3dsXpi5fw8No9OHjWIHVbRD0CwxERUScV5O2KzXNiERHggcp6E2a8k409xZVSt0XU7TEcERF1Yr7uGqTOHodxoV6oN5ox6929+OpAmdRtEXVrDEdERJ2ch1aF9389BpOH+MNkseIPH/2Af+eclLotom6L4YiIqAvQqhRY9XgkHh0TBCGARekH8ff/FvFu2kROwHBERNRFKOQyvPLgUMy9ayAA4M2MY3jps0OwWhmQiNoTwxERURcik8nwbFwYXro/AjIZ8EHWSTyzIR8mM39uhKi9MBwREXVBs34WgpUzRkGlkOHzgnP47Qd70WA0S90WUbfAcERE1EU9MCIQ62eOhqtagZ1FlXhsXTaq6o1St0XU5TEcERF1YXfc5ouPnhiH3q4qFJwx4JG3s3CmulHqtoi6NIYjIqIubmR/T2yaE4u+ni4oudCAh9dk4dj5OqnbIuqyGI6IiLqBgX16YfOTMRjUpxfKay/jkbVZ2HfyotRtEXVJDEdERN1EgM4Fm+bEIDLIE4ZLTXj8Hzn49sh5qdsi6nIYjoiIuhFPVzX+/btx+HmYLy43WfHEh/uQtu+M1G0RdSkMR0RE3YyLWoF3fhWN6aP6wmIVeG5TAdZ9VyJ1W0RdBsMREVE3pFLI8fojI/C720MAAH/7z2Es/eowf26EqA0YjoiIuim5XIZFvxiM5CnhAIC3M0vw4ub9MFt4N22i1jAcERF1YzKZDHMmDMCrDw+HXAZs2ncGc/61D5ebLFK3RtRpMRwREfUA8dH98XZiNDRKOb45XIHE9TkwNDZJ3RZRp8RwRETUQ9wT4Yd//nYs3LVK7D1RjYR3snC+9rLUbRF1OgxHREQ9yJgQL2ycHYM+7hocKa/DQ2v2oORCvdRtEXUqDEdERD3M4AAPpD0Zi2BvV5ypvoRH1mbhwBmD1G0RdRoMR0REPVB/L1dsfjIWQ/t6oKrBhBnvZGH38Uqp2yLqFBiOiIh6KJ9eGnz8xDjEDvBGg8mCX7+3F1/uL5O6LSLJMRwREfVg7loV3vv1aNw7zB8mixVPffwD/pl9Uuq2iCTFcERE1MNplAr8v0cj8fjYIAgB/GnLQaz45hjvpk09FsMRERFBIZfhr9OG4pm7BwEAVnxThD9/eggWKwMS9TwMR0REBODK3bTn3XMb/m/qEMhkwD+zT2Juah6MZt5Nm3oWhiMiIrLzq5hg/L9HR0GlkOHL/WX4zft7UW80S90WUYdhOCIiombuGx6I92aNgatagd3Hq/DoO9moqjdK3RZRh2A4IiKiFt0+yAepvx8HLzc1Dpw14OG1WTh9sVHqtoicjuGIiIiua3g/T2yeE4O+ni4orWzAQ2v24Eh5rdRtETnVTYWj1atXIyQkBFqtFlFRUdi5c2er9ZmZmYiKioJWq0VoaCjWrl3brCYtLQ0RERHQaDSIiIhAenq6w9s9f/48Zs2ahcDAQLi6umLy5MkoKipqNk9WVhbuuusuuLm5wdPTE3feeScuXbrk4LNARNQzhPr2QtqTsQjzc0dFnRHxa7Ow98RFqdsichqHw9GGDRuQlJSERYsWIS8vD+PHj8eUKVNw6tSpFutLS0tx7733Yvz48cjLy8PChQsxd+5cpKWl2WqysrKQkJCAxMREFBQUIDExEfHx8cjJyWnzdoUQmDZtGkpKSvDpp58iLy8Per0eEydORENDg922Jk+ejLi4OHz//ffYu3cvnnrqKcjlPIlGRHQ9/jotNs6OQbS+N2ovm/HLf+Tgm8LzUrdF5BQy4eBdvsaOHYvIyEisWbPGNjZ48GBMmzYNS5cubVY/f/58fPbZZzh8+LBtbM6cOSgoKEBWVhYAICEhAbW1tfjqq69sNZMnT0bv3r3x8ccft2m7x44dQ1hYGA4ePIghQ4YAACwWC/r06YNly5bhd7/7HQBg3LhxuOeee/Dyyy87sts2tbW10Ol0MBgM8PDwuKk5iIi6qksmC5766Af890gFFHIZUqYPwyPR/aVui+iGHHn/duh0iclkwr59+xAXF2c3HhcXhz179rS4TlZWVrP6SZMmITc3F01NTa3WXJ2zLds1Gq98i0Kr1dqWKxQKqNVq7Nq1CwBQUVGBnJwc9OnTB7GxsfDz88OECRNsy1tiNBpRW1tr9yAi6qlc1AqsTYzCQ5H9YLEKvLB5P97OLJa6LaJ25VA4qqyshMVigZ+fn924n58fysvLW1ynvLy8xXqz2YzKyspWa67O2ZbthoeHQ6/XY8GCBaiurobJZEJKSgrKy8tRVnblhxRLSkoAAC+99BKeeOIJbN26FZGRkbj77rtbvDYJAJYuXQqdTmd79O/P/yERUc+mUsjx+iPDMfuOUADA0q+O4JX/HIaVd9OmbuKmLrSRyWR2fwshmo3dqP6n422Zs7UalUqFtLQ0HDt2DF5eXnB1dcWOHTswZcoUKBQKAIDVagUAzJ49G7/+9a8xatQoLF++HGFhYXj33Xdb7H3BggUwGAy2x+nTp6+7n0REPYVMJsOCewdjwZRwAMA735Xghc370WSxStwZ0a1TOlLs4+MDhULR7CxRRUVFs7M6V/n7+7dYr1Qq4e3t3WrN1Tnbut2oqCjk5+fDYDDAZDLB19cXY8eORXR0NAAgICAAABAREWE3z+DBg697QblGo4FGo2lxGRFRTzd7wgB499Jgftp+pP1wBtWNJqx6LBIuaoXUrRHdNIfOHKnVakRFRSEjI8NuPCMjA7GxsS2uExMT06x+27ZtiI6OhkqlarXm6pyOblen08HX1xdFRUXIzc3F1KlTAQDBwcEIDAzE0aNH7eqPHTsGvV5/o90nIqIWPBzVD+8kRkGjlOPbIxX45focGBqbpG6L6OYJB6WmpgqVSiXWr18vCgsLRVJSknBzcxMnTpwQQgiRnJwsEhMTbfUlJSXC1dVVzJs3TxQWFor169cLlUolNm/ebKvZvXu3UCgUIiUlRRw+fFikpKQIpVIpsrOz27xdIYTYuHGj2L59uyguLhZbtmwRer1eTJ8+3a7/5cuXCw8PD7Fp0yZRVFQkFi9eLLRarTh+/Hib9t9gMAgAwmAwOPrUERF1a3tLq8SwJVuFfv4X4p43d4iymktSt0Rk48j7t8PhSAghVq1aJfR6vVCr1SIyMlJkZmbals2cOVNMmDDBrn7Hjh1i1KhRQq1Wi+DgYLFmzZpmc27atEmEhYUJlUolwsPDRVpamkPbFUKIlStXin79+gmVSiWCgoLE4sWLhdFobDbP0qVLRb9+/YSrq6uIiYkRO3fubPO+MxwREV3fkbJaMeZvGUI//wsRu/S/4nhFndQtEQkhHHv/dvg+Rz0d73NERNS60xcbMfPd71FS2QAvNzXemzUaI/p7St0W9XBOu88RERHRjfT3csWmOTEY3k+Hiw0mPLouGzuLLkjdFlGbMRwREVG78+6lwUdPjMPtA33QaLLgN+/vxecF56Rui6hNGI6IiMgpemmUWD8rGr8YHoAmi8Dc1Dx8mHVC6raIbojhiIiInEajVODvM0bhVzF6CAH8+dNDeDPjGHi5K3VmDEdERORUCrkMf3lgCOZNvA0A8Pf/FmHxloOw8OdGqJNiOCIiIqeTyWR4ZuIg/HXaUMhkwL9zTuHpj3+A0WyRujWiZhiOiIiow/xynB6rHouEWiHHfw6U49fv7UXdZd5NmzoXhiMiIupQ9w4LwPu/Hg03tQJ7iqvw6LpsVNYbpW6LyIbhiIiIOlzsQB+k/j4G3m5qHDxbi4fX7MHpi41St0UEgOGIiIgkMqyfDpufjEW/3i44UdWI6Wv24HBZrdRtETEcERGRdEJ83JD2ZCzC/d1xoc6I+Lez8H3pRanboh6O4YiIiCTl56HFhtkxGB3cG3WXzUhcn4OMwvNSt0U9GMMRERFJTueiwj9/OxYTB/eB0WzFnH/tw8bc01K3RT0UwxEREXUKWpUCa38ZhUei+sFiFXhx836szSzm3bSpwzEcERFRp6FUyPHqw8MxZ8IAAEDKV0fwty8Pw8q7aVMHYjgiIqJORSaTIXlKOBb/YjAA4B+7SvH8pgI0WawSd0Y9hVLqBoiIiFryu/Gh8HJT44XN+/FJ3llUN5rwp/siUF57GSE+bgjQuUjdInVTDEdERNRpTY/sh96uajz5733YfvQCth/NBADIZcDS6cOQMDpI4g6pO+LHakRE1Kn9PLwP/j5jlN2YVQALPzmIMsMlibqi7ozhiIiIOr1e2uYfdFiEwIlK/uQItT+GIyIi6vRCfNwgl9mPyWRAsI+rNA1Rt8ZwREREnV6AzgVLpw+DQvZjQlIr5OA3/MkZGI6IiKhLSBgdhF3JP8e/fjsGEQEeMJqteG5jPu+BRO2O4YiIiLqMAJ0Lbh/ki1WPR8JVrUB2yUWs21kidVvUzTAcERFRlxPi44Y/3xcBAHh921EcOmeQuCPqThiOiIioS0oY3R9xEX5osggkpebjcpNF6paom2A4IiKiLkkmkyHloeHwddegqKIeKV8dkbol6iYYjoiIqMvyclPjtYeHAwDe33MCmccuSNwRdQcMR0RE1KXdGdYHs2KDAQDPbyrAxQaTtA1Rl8dwREREXV7ylHAM7NMLF+qMSE7bDyH49X66eQxHRETU5WlVCqycMRIqhQzbCs9jY+5pqVuiLozhiIiIuoUhgTo8HxcGAPjL54U4UdkgcUfUVTEcERFRt/G78aEYF+qFRpMFSRvyYbZYpW6JuiCGIyIi6jYUchneiB8Jd60S+adr8P++PS51S9QFMRwREVG30tfTBX97cBgA4K3tx7HvZLXEHVFXw3BERETdzgMjAjFtZCAsVoFnN+aj3miWuiXqQhiOiIioW/rL1KHo6+mCk1WN+L/PD0ndDnUhDEdERNQt6VxUeDN+BGQyYGPuGWw9WCZ1S9RFMBwREVG3NTbUG3MmDAAAJH9yAOdrL0vcEXUFDEdERNStzZt4G4b29UBNYxOe31QAq5V3z6bWMRwREVG3plbKsSJhFLQqOXYWVeKDrBNSt0SdHMMRERF1ewP79MKiewcDAJZ+dQTHztdJ3BF1ZgxHRETUI/xynB4/D/OFyWzF3I/zYDRbpG6JOimGIyIi6hFkMhlefXgEvN3UOFJehze2HZO6JeqkGI6IiKjH8HXXIOWh4QCAdTtLsOd4pcQdUWfEcERERD3KPRF+eHRMEIQAnttUAENjk9QtUSfDcERERD3On+4bjBAfN5QZLmPhlgMQgl/vpx8xHBERUY/jqlZiRcJIKOQyfLm/DFvyz0rdEnUiDEdERNQjjejviaS7BwEA/rzlEE5fbJS4I+osGI6IiKjH+sPPByJa3xt1RjOe3ZgPC++eTWA4IiKiHkwhl2F5wkj00iix90Q11mYWS90SdQIMR0RE1KP193LFSw8MAQAszziG/WdqpG2IJMdwREREPd5DkX1x7zB/mK0CSan5aDSZpW6JJMRwREREPZ5MJsMrDw6Dn4cGJZUNeOU/h6VuiSTEcERERATA01WNNx4ZCQD4V/Yp/PfweWkbIskwHBEREf3P7YN88NvbQwAAL27ejwt1Rok7IikwHBEREV3jhUlhCPd3R1WDCclp+3n37B6I4YiIiOgaWpUCK2aMhFohx3+PVODfOaekbok62E2Fo9WrVyMkJARarRZRUVHYuXNnq/WZmZmIioqCVqtFaGgo1q5d26wmLS0NERER0Gg0iIiIQHp6usPbPX/+PGbNmoXAwEC4urpi8uTJKCoqarEnIQSmTJkCmUyGLVu2tH3niYio2wv398CLk8MAAH/9shDFF+ol7og6ksPhaMOGDUhKSsKiRYuQl5eH8ePHY8qUKTh1quVkXVpainvvvRfjx49HXl4eFi5ciLlz5yItLc1Wk5WVhYSEBCQmJqKgoACJiYmIj49HTk5Om7crhMC0adNQUlKCTz/9FHl5edDr9Zg4cSIaGhqa9bVixQrIZDJHd5+IiHqI3/wsBLcP9MHlJiuSUvPRZLFK3RJ1FOGgMWPGiDlz5tiNhYeHi+Tk5BbrX3zxRREeHm43Nnv2bDFu3Djb3/Hx8WLy5Ml2NZMmTRIzZsxo83aPHj0qAIiDBw/alpvNZuHl5SXWrVtnt15+fr7o16+fKCsrEwBEenr6Dfb6RwaDQQAQBoOhzesQEVHXVFZzSQx/6Wuhn/+FeHXrYanboVvgyPu3Q2eOTCYT9u3bh7i4OLvxuLg47Nmzp8V1srKymtVPmjQJubm5aGpqarXm6pxt2a7ReOUbBVqt1rZcoVBArVZj165dtrHGxkY8+uijeOutt+Dv73/DfTYajaitrbV7EBFRz+Cv02Lp9GEAgNU7ivF96UWJO6KO4FA4qqyshMVigZ+fn924n58fysvLW1ynvLy8xXqz2YzKyspWa67O2ZbthoeHQ6/XY8GCBaiurobJZEJKSgrKy8tRVlZmW2fevHmIjY3F1KlT27TPS5cuhU6nsz369+/fpvWIiKh7uHdYAB6O6gchgHkb8lF7uUnqlsjJbuqC7J9eqyOEaPX6nZbqfzreljlbq1GpVEhLS8OxY8fg5eUFV1dX7NixA1OmTIFCoQAAfPbZZ/j222+xYsWKNuzlFQsWLIDBYLA9Tp8+3eZ1iYioe1hyfwT6e7ngbM0lvPTpIanbISdzKBz5+PhAoVA0O0tUUVHR7KzOVf7+/i3WK5VKeHt7t1pzdc62bjcqKgr5+fmoqalBWVkZtm7diqqqKoSEXLmh17fffovi4mJ4enpCqVRCqVQCAB566CHceeedLfav0Wjg4eFh9yAiop7FXavC8viRkMuAT/LO4vOCc1K3RE7kUDhSq9WIiopCRkaG3XhGRgZiY2NbXCcmJqZZ/bZt2xAdHQ2VStVqzdU5Hd2uTqeDr68vioqKkJuba/sILTk5Gfv370d+fr7tAQDLly/He++918ZngYiIeqLoYC889fOBAIBF6QdQZrgkcUfkNI5e7Z2amipUKpVYv369KCwsFElJScLNzU2cOHFCCCFEcnKySExMtNWXlJQIV1dXMW/ePFFYWCjWr18vVCqV2Lx5s61m9+7dQqFQiJSUFHH48GGRkpIilEqlyM7ObvN2hRBi48aNYvv27aK4uFhs2bJF6PV6MX369Fb3B/y2GhERtZHJbBEP/L+dQj//C/HoO1nCYrFK3RK1kSPv3w6HIyGEWLVqldDr9UKtVovIyEiRmZlpWzZz5kwxYcIEu/odO3aIUaNGCbVaLYKDg8WaNWuazblp0yYRFhYmVCqVCA8PF2lpaQ5tVwghVq5cKfr16ydUKpUICgoSixcvFkajsdV9YTgiIiJHFFfUifDFXwn9/C/EO5nFUrdDbeTI+7dMCP5ojCNqa2uh0+lgMBh4/RERUQ/18fensOCTA1Ar5Njyx58hIpDvB52dI+/f/G01IiIiB80Y3R8TB/vBZLEiaUMeLjdZpG6J2hHDERERkYNkMhmWPTQMPr00OHa+Hsu2HpG6JWpHDEdEREQ3wbuXBq89MhwA8N7uE/ju2AWJO6L2wnBERER0k34e1ge/itEDAJ7fVICLDSaJO6L2wHBERER0CxZMGYwBvm6oqDNi4ScHwO85dX0MR0RERLfARa3AyhmjoFLIsPVQOTbtOyN1S3SLGI6IiIhu0dC+Ojx7TxgA4C+fHcLJqgaJO6JbwXBERETUDn5/RyjGhHihwWRB0oZ8mC1WqVuim8RwRERE1A4UchmWJ4yEu1aJvFM1WLW9WOqW6CYxHBEREbWTvp4u+Ou0oQCAv39bhLxT1RJ3RDeD4YiIiKgdTR3ZFw+MCITFKjBvQz4ajGapWyIHMRwRERG1s5enDUWgTosTVY14+YtCqdshBzEcERERtTOdiwpvxI+ETAak7j2Nrw+VS90SOYDhiIiIyAliBnjj93eEAgCS0/ajovayxB1RWzEcEREROclz94QhIsAD1Y1NeH7zft49u4tgOCIiInIStVKOlTNGQqOU47tjF/Bh1kmpW6I2YDgiIiJyokF+7lh472AAwCv/OYyi83USd0Q3wnBERETkZL+K0WPCbb4wmq14JjUfRrNF6paoFQxHRERETiaTyfDaI8Ph5aZGYVkt3sw4JnVL1AqGIyIiog7Qx12LlOnDAADvfFeCrOIqiTui62E4IiIi6iBxQ/zx6Jj+EAJ4dmM+DI1NUrdELWA4IiIi6kCLfxGBYG9XlBku40+fHpS6HWoBwxEREVEHctMosTxhJBRyGT4rOIcteWelbol+guGIiIiog40K6o1n7h4EAPjTloM4U90ocUd0LYYjIiIiCfzhzgGIDPJEndGMZzcWwGLl3bM7C4YjIiIiCSgVcqxIGAU3tQLfl17E298VS90S/Q/DERERkUSCvF3x0gNDAABvbjuGg2cNEndEAMMRERGRpB6O6ocpQ/1htgo8k5qHSybePVtqDEdEREQSkslkeOXBYejjrkHxhQa88p/DUrfU4zEcERERSay3mxpvxI8AAPwz+yS2H6mQuKOejeGIiIioExg/yBe/+VkIAOCFzQWorDdK3FHPxXBERETUSbw4OQy3+fVCZb0JyWn7IQS/3i8FhiMiIqJOQqtSYOWMUVAr5PjmcAU+/v601C31SAxHREREncjgAA+8ODkMAPDyF4UouVAvcUc9D8MRERFRJ/Obn4UgdoA3LjVZMG9DPposVqlb6lEYjoiIiDoZuVyGN+JHQOeiQsEZA/7+3yKpW+pRGI6IiIg6oQCdC155cBgAYNX248g9cVHijnoOhiMiIqJO6hfDAzA9si+sAkjakI+6y01St9QjMBwRERF1Yn95YAj69XbBmepLeOmzQqnb6REYjoiIiDoxd60KyxNGQi4D0n44gy/3l0ndUrfHcERERNTJjQ72wh/uHAgAWJh+AGWGSxJ31L0xHBEREXUBz0wchOH9dDBcasLzmwpgtfLu2c7CcERERNQFqBRyrEgYCReVAruPV+Hd3aVSt9RtMRwRERF1EaG+vbD4vsEAgFe3HsXhslqJO+qeGI6IiIi6kMfGBGHi4D4wWaxISs3H5SaL1C11OwxHREREXYhMJkPKQ8Ph00uNo+fr8NrXR6VuqdthOCIiIupifHpp8OrDwwEA63eVYmfRBYk76l4YjoiIiLqgu8L98MtxQQCA5zcVoLrBJHFH3QfDERERURe16N4IhPq64XytEQvTD0AIfr2/PTAcERERdVEuagVWJoyCUi7DVwfLsXnfGalb6hYYjoiIiLqwYf10mHfPbQCAlz47hFNVjRJ31PUxHBEREXVxcyYMwJhgLzSYLJi3MR9mi1Xqlro0hiMiIqIuTiGX4Y34EXDXKLHvZDVW7yiWuqUujeGIiIioG+jv5Yr/mzYEALDyv0XIP10jbUNdGMMRERFRNzFtZF/cPyIQFqtAUmoeGoxmqVvqkhiOiIiIugmZTIa/Th2KAJ0WJ6oa8dcvC6VuqUtiOCIiIupGdK4qvBE/AjIZ8PH3p7HtULnULXU5DEdERETdTOwAH/x+fCgAIPmTA6iouyxxR13LTYWj1atXIyQkBFqtFlFRUdi5c2er9ZmZmYiKioJWq0VoaCjWrl3brCYtLQ0RERHQaDSIiIhAenq6w9s9f/48Zs2ahcDAQLi6umLy5MkoKiqyLb948SKefvpphIWFwdXVFUFBQZg7dy4MBsPNPA1ERESd1rNxt2FwgAcuNpjwwqb9vHu2AxwORxs2bEBSUhIWLVqEvLw8jB8/HlOmTMGpU6darC8tLcW9996L8ePHIy8vDwsXLsTcuXORlpZmq8nKykJCQgISExNRUFCAxMRExMfHIycnp83bFUJg2rRpKCkpwaeffoq8vDzo9XpMnDgRDQ0NAIBz587h3LlzeP3113HgwAG8//772Lp1K3772986+jQQERF1ahqlAitnjIRGKUfmsQv4Z/ZJqVvqOoSDxowZI+bMmWM3Fh4eLpKTk1usf/HFF0V4eLjd2OzZs8W4ceNsf8fHx4vJkyfb1UyaNEnMmDGjzds9evSoACAOHjxoW242m4WXl5dYt27ddfdn48aNQq1Wi6ampuvWXMtgMAgAwmAwtKmeiIhISu/tKhH6+V+I2xb9RxSdr5W6Hck48v7t0Jkjk8mEffv2IS4uzm48Li4Oe/bsaXGdrKysZvWTJk1Cbm4umpqaWq25Omdbtms0GgEAWq3WtlyhUECtVmPXrl3X3SeDwQAPDw8olcoWlxuNRtTW1to9iIiIuopfxQTjjtt8YTRb8UxqPkxm3j37RhwKR5WVlbBYLPDz87Mb9/PzQ3l5y1fDl5eXt1hvNptRWVnZas3VOduy3fDwcOj1eixYsADV1dUwmUxISUlBeXk5ysrKWuytqqoKL7/8MmbPnn3dfV66dCl0Op3t0b9//+vWEhERdTZyuQyvPzwcvV1VOHSuFm9mHJO6pU7vpi7Ilslkdn8LIZqN3aj+p+NtmbO1GpVKhbS0NBw7dgxeXl5wdXXFjh07MGXKFCgUimY91dbW4he/+AUiIiKwZMmS6/a+YMECGAwG2+P06dPXrSUiIuqM+nhosXT6cADA298VI7ukSuKOOjeHwpGPjw8UCkWzs0QVFRXNzupc5e/v32K9UqmEt7d3qzVX52zrdqOiopCfn4+amhqUlZVh69atqKqqQkhIiN16dXV1mDx5Mnr16oX09HSoVKrr7rNGo4GHh4fdg4iIqKuZPNQfCdH9IQTw3MYCGC41Sd1Sp+VQOFKr1YiKikJGRobdeEZGBmJjY1tcJyYmpln9tm3bEB0dbQsl16u5Oqej29XpdPD19UVRURFyc3MxdepU27La2lrExcVBrVbjs88+s7tGiYiIqDv78/0R0Hu74mzNJfz504NSt9N5OXq1d2pqqlCpVGL9+vWisLBQJCUlCTc3N3HixAkhhBDJyckiMTHRVl9SUiJcXV3FvHnzRGFhoVi/fr1QqVRi8+bNtprdu3cLhUIhUlJSxOHDh0VKSopQKpUiOzu7zdsV4so3z7Zv3y6Ki4vFli1bhF6vF9OnT7ctr62tFWPHjhXDhg0Tx48fF2VlZbaH2Wxu0/7z22pERNSV7Tt5UYQu+FLo538htuSdkbqdDuPI+7fD4UgIIVatWiX0er1Qq9UiMjJSZGZm2pbNnDlTTJgwwa5+x44dYtSoUUKtVovg4GCxZs2aZnNu2rRJhIWFCZVKJcLDw0VaWppD2xVCiJUrV4p+/foJlUolgoKCxOLFi4XRaLQt3759uwDQ4qO0tLRN+85wREREXd3yjKNCP/8LMXTJVnGmulHqdjqEI+/fMiF4y0xH1NbWQqfT2W4BQERE1NWYLVY88nYW8k7VYGyIFz56YhwU8ut/sao7cOT9m7+tRkRE1MMoFXKsSBgJV7UCOaUXsW5nidQtdSoMR0RERD2Q3tsNL90/BADwxrajOHiWvzN6FcMRERFRD/VIdD9MGuKHJotA0oZ8XDJZpG6pU2A4IiIi6qFkMhmWTh+OPu4aHK+oR8pXh6VuqVNgOCIiIurBvNzUeO2REQCAD7JOYvvRCok7kh7DERERUQ834TZfzIoNBgC8sGk/quqN0jYkMYYjIiIiQvKUcNzm1wuV9UYkf3IAPflOPwxHREREBK1KgRUJo6BWyJFReB4b9vbcH1pnOCIiIiIAQESgB56fdBsA4C+fF6K0skHijqTBcEREREQ2v7s9FDGh3rjUZEHShnw0WaxSt9ThGI6IiIjIRi6X4Y34EfDQKlFwugb/79vjUrfU4RiOiIiIyE6gpwv+9uAwAMBb3xZh38mLEnfUsRiOiIiIqJn7RwTiwVF9YRVA0oZ81BvNUrfUYRiOiIiIqEV/mToEfT1dcPriJfzls0NSt9NhGI6IiIioRR5aFZYnjIRMBmzadwZfHSiTuqUOwXBERERE1zUmxAtPThgAAFiQfgDlhssSd+R8DEdERETUqqSJt2FYXx1qGpvwwuYCWK3d++7ZDEdERETUKrVSjuUJI6FVybGzqBLv7TkhdUtOxXBERERENzSwTy8s+kUEAGDZ1iM4Ul4rcUfOw3BEREREbfLLsUG4K7wPTGYrklLzcbnJInVLTsFwRERERG0ik8mw7KHh8HZT40h5HV7/+qjULTkFwxERERG1ma+7Bq8+PBwA8I9dpdh9vFLijtofwxERERE55O7Bfnh8bBAA4LmNBahpNEncUftiOCIiIiKHLfrFYIT6uKG89jIWpR+EEN3n6/0MR0REROQwV7USK2aMhFIuw5cHyvDJD2elbqndMBwRERHRTRnezxPz7rkNALDks0M4fbFR4o7aB8MRERER3bQ5EwYgWt8b9UYz5m3Ih9lilbqlW8ZwRERERDdNIZdhecJI9NIokXuyGmszi6Vu6ZYxHBEREdEt6e/liv+bOgQAsOKbIhScrpG2oVvEcERERES37MFRffGL4QEwWwXmbchHo8ksdUs3jeGIiIiIbplMJsPfpg2Fv4cWJZUN+OuXh6Vu6aYxHBEREVG78HRV4834EQCAj3JO4ZvC8xJ3dHMYjoiIiKjdxA70wRPjQwAA89P240KdUeKOHMdwRERERO3q+UlhCPd3R1WDCS9uLuhyd89mOCIiIqJ2pVEqsHLGKKiVcmw/egH/yjkldUsOYTgiIiKidhfm747kyeEAgL99WYjjFfUSd9R2DEdERETkFLNigzF+kA8uN1mRtCEPJnPXuHs2wxERERE5hVwuw+uPjICnqwoHz9ZixTfHpG6pTRiOiIiIyGn8PLRImT4MALAmsxg5JVUSd3RjDEdERETkVJOHBuCRqH4QAnh2YwFqLzdJ3VKrGI6IiIjI6ZY8MARBXq44W3MJSz49JHU7rWI4IiIiIqfrpVFiecJIyGVAet5ZfFZwTuqWrovhiIiIiDpElL43nrprEABgUfoBnK25JHFHLWM4IiIiog7z9F0DMbK/J+oum/HcxnxYrZ3v7tkMR0RERNRhVAo5lieMhKtageySi/jHrhKpW2qG4YiIiIg6VIiPG/58XwQA4LWvj+LQOYPEHdljOCIiIqIOlzC6P+Ii/NBkEUhKzcflJovULdkwHBEREVGHk8lkSHloOHzdNSiqqEfKV0ekbsmG4YiIiIgk4eWmxmsPDwcAvL/nBDKPXZC4oysYjoiIiEgyd4b1wazYYADA85sKcLisFnuKK1FmkO5r/krJtkxEREQEIHlKOHYdr8TxinpMWbkTACCXAUunD0PC6KAO74dnjoiIiEhSWpUCf/rFYLsxqwAWfnJQkjNIDEdEREQkOZWyeSSxCIETlY0d3gvDEREREUkuxMcNcpn9mEImQ7CPa4f3wnBEREREkgvQuWDp9GFQyK4kJIVMhlemD0WAzqXDe+EF2URERNQpJIwOwh23+eJEZSOCfVwlCUYAwxERERF1IgE6F8lC0VX8WI2IiIjoGjcVjlavXo2QkBBotVpERUVh586drdZnZmYiKioKWq0WoaGhWLt2bbOatLQ0REREQKPRICIiAunp6Q5v9/z585g1axYCAwPh6uqKyZMno6ioyK7GaDTi6aefho+PD9zc3PDAAw/gzJkzN/EsEBERUXfkcDjasGEDkpKSsGjRIuTl5WH8+PGYMmUKTp061WJ9aWkp7r33XowfPx55eXlYuHAh5s6di7S0NFtNVlYWEhISkJiYiIKCAiQmJiI+Ph45OTlt3q4QAtOmTUNJSQk+/fRT5OXlQa/XY+LEiWhoaLDNk5SUhPT0dKSmpmLXrl2or6/HfffdB4ul8/zgHREREUlIOGjMmDFizpw5dmPh4eEiOTm5xfoXX3xRhIeH243Nnj1bjBs3zvZ3fHy8mDx5sl3NpEmTxIwZM9q83aNHjwoA4uDBg7blZrNZeHl5iXXr1gkhhKipqREqlUqkpqbaas6ePSvkcrnYunXrDfddCCEMBoMAIAwGQ5vqiYiISHqOvH87dObIZDJh3759iIuLsxuPi4vDnj17WlwnKyurWf2kSZOQm5uLpqamVmuuztmW7RqNRgCAVqu1LVcoFFCr1di1axcAYN++fWhqarKbJzAwEEOHDr1u/0ajEbW1tXYPIiIi6r4cCkeVlZWwWCzw8/OzG/fz80N5eXmL65SXl7dYbzabUVlZ2WrN1Tnbst3w8HDo9XosWLAA1dXVMJlMSElJQXl5OcrKymzbUavV6N27d5v7X7p0KXQ6ne3Rv3//6z4/RERE1PXd1AXZMpn9LSyFEM3GblT/0/G2zNlajUqlQlpaGo4dOwYvLy+4urpix44dmDJlChQKRav701r/CxYsgMFgsD1Onz7d6lxERETUtTl0nyMfHx8oFIpmZ1kqKiqandW5yt/fv8V6pVIJb2/vVmuuztnW7UZFRSE/Px8GgwEmkwm+vr4YO3YsoqOjbdsxmUyorq62O3tUUVGB2NjYFvvXaDTQaDTXfU6IiIioe3HozJFarUZUVBQyMjLsxjMyMq4bLmJiYprVb9u2DdHR0VCpVK3WXJ3T0e3qdDr4+vqiqKgIubm5mDp1KoAr4UmlUtnNU1ZWhoMHD163fyIiIuphHL3aOzU1VahUKrF+/XpRWFgokpKShJubmzhx4oQQQojk5GSRmJhoqy8pKRGurq5i3rx5orCwUKxfv16oVCqxefNmW83u3buFQqEQKSkp4vDhwyIlJUUolUqRnZ3d5u0KIcTGjRvF9u3bRXFxsdiyZYvQ6/Vi+vTpdv3PmTNH9OvXT3zzzTfihx9+EHfddZcYMWKEMJvNbdp/fluNiIio63Hk/dvhcCSEEKtWrRJ6vV6o1WoRGRkpMjMzbctmzpwpJkyYYFe/Y8cOMWrUKKFWq0VwcLBYs2ZNszk3bdokwsLChEqlEuHh4SItLc2h7QohxMqVK0W/fv2ESqUSQUFBYvHixcJoNNrVXLp0STz11FPCy8tLuLi4iPvuu0+cOnWqzfvOcERERNT1OPL+LRPif1dHU5sYDAZ4enri9OnT8PDwkLodIiIiaoPa2lr0798fNTU10Ol0rdbyh2cdVFdXBwD8Sj8REVEXVFdXd8NwxDNHDrJarTh37hzc3d1bvX3BzbiaanlWqnPhcel8eEw6Jx6XzonH5QohBOrq6hAYGAi5vPXvo/HMkYPkcjn69evn1G14eHj06H/AnRWPS+fDY9I58bh0TjwuuOEZo6tu6iaQRERERN0VwxERERHRNRiOOhGNRoMlS5bwjtydDI9L58Nj0jnxuHROPC6O4wXZRERERNfgmSMiIiKiazAcEREREV2D4YiIiIjoGgxHRERERNdgOGondXV1SEpKgl6vh4uLC2JjY7F3797r1u/YsQMymazZ48iRI3Z1aWlpiIiIgEajQUREBNLT05vNtXr1aoSEhECr1SIqKgo7d+5s9/3rqqQ6Li+99FKzOfz9/Z2yj12NM47JoUOH8NBDDyE4OBgymQwrVqxocS6+Vq5PquPC10rrnHFc1q1bh/Hjx6N3797o3bs3Jk6ciO+//77ZXD359cJw1E5+97vfISMjA//85z9x4MABxMXFYeLEiTh79myr6x09ehRlZWW2x6BBg2zLsrKykJCQgMTERBQUFCAxMRHx8fHIycmx1WzYsAFJSUlYtGgR8vLyMH78eEyZMgWnTp1y2r52JVIdFwAYMmSI3RwHDhxwyj52Nc44Jo2NjQgNDUVKSsp131j5WmmdVMcF4GulNc44Ljt27MCjjz6K7du3IysrC0FBQYiLi7Obs8e/XgTdssbGRqFQKMQXX3xhNz5ixAixaNGiFtfZvn27ACCqq6uvO298fLyYPHmy3dikSZPEjBkzbH+PGTNGzJkzx64mPDxcJCcnO7gX3Y+Ux2XJkiVixIgRN917d+WsY3ItvV4vli9f3mycr5Xrk/K48LVyfR1xXIQQwmw2C3d3d/HBBx/Yxnr664VnjtqB2WyGxWKBVqu1G3dxccGuXbtaXXfUqFEICAjA3Xffje3bt9sty8rKQlxcnN3YpEmTsGfPHgCAyWTCvn37mtXExcXZanoyqY7LVUVFRQgMDERISAhmzJiBkpKSW9ib7sFZx+RG+FppnVTH5Sq+VlrWUcelsbERTU1N8PLyAsDXC8CP1dqFu7s7YmJi8PLLL+PcuXOwWCz417/+hZycHJSVlbW4TkBAAN555x2kpaXhk08+QVhYGO6++2589913tpry8nL4+fnZrefn54fy8nIAQGVlJSwWS6s1PZlUxwUAxo4diw8//BBff/011q1bh/LycsTGxqKqqso5O9tFOOuY3AhfK62T6rgAfK20pqOOS3JyMvr27YuJEycC4OsFAD9Way/Hjx8Xd9xxhwAgFAqFGD16tHj88cfF4MGD2zzHfffdJ+6//37b3yqVSnz00Ud2Nf/617+ERqMRQghx9uxZAUDs2bPHruavf/2rCAsLu4W96T6kOC4tqa+vF35+fuKNN95wfCe6GWcck2u19PENXys3JsVxaQlfK/acfVyWLVsmevfuLQoKCmxjfL3wY7V2M2DAAGRmZqK+vh6nT5/G999/j6amJoSEhLR5jnHjxqGoqMj2t7+/f7OUXlFRYUvzPj4+UCgUrdb0dFIcl5a4ublh2LBhdvP0VM44JjfC18qNSXFcWsLXij1nHpfXX38dr7zyCrZt24bhw4fbxvl64cdq7c7NzQ0BAQGorq7G119/jalTp7Z53by8PAQEBNj+jomJQUZGhl3Ntm3bEBsbCwBQq9WIiopqVpORkWGroSs68ri0xGg04vDhw3bz9HTteUxuhK+VtuvI49ISvlZa1t7H5bXXXsPLL7+MrVu3Ijo62m4ZXy/gx2rtZevWreKrr74SJSUlYtu2bWLEiBFizJgxwmQyCSGESE5OFomJibb65cuXi/T0dHHs2DFx8OBBkZycLACItLQ0W83u3buFQqEQKSkp4vDhwyIlJUUolUqRnZ1tq0lNTRUqlUqsX79eFBYWiqSkJOHm5iZOnDjRcTvfiUl1XJ577jmxY8cOUVJSIrKzs8V9990n3N3deVyEc46J0WgUeXl5Ii8vTwQEBIjnn39e5OXliaKiIlsNXyutk+q48LXSOmccl2XLlgm1Wi02b94sysrKbI+6ujpbTU9/vTActZMNGzaI0NBQoVarhb+/v/jjH/8oampqbMtnzpwpJkyYYPt72bJlYsCAAUKr1YrevXuL22+/XXz55ZfN5t20aZMICwsTKpVKhIeH2/0Dv2rVqlVCr9cLtVotIiMjRWZmplP2sSuS6rgkJCSIgIAAoVKpRGBgoJg+fbo4dOiQ0/azK3HGMSktLRUAmj2unUcIvlZaI9Vx4Wuldc44Lnq9vsXjsmTJEru6nvx6kQkhRMeeqyIiIiLqvHjNEREREdE1GI6IiIiIrsFwRERERHQNhiMiIiKiazAcEREREV2D4YiIiIjoGgxHRERERNdgOCIiIiK6BsMRERER0TUYjoiIiIiuwXBEREREdA2GIyIiIqJr/H+fdaIbCH6vEwAAAABJRU5ErkJggg==", "text/plain": [ "
" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "ds = xr.open_zarr(\"SwashParticles.zarr\")\n", "\n", "plt.plot(ds.lon.T, ds.lat.T, \".-\")\n", "plt.show()" ] }, { "attachments": {}, "cell_type": "markdown", "metadata": {}, "source": [ "Note that, even though we use 2-dimensional `AdvectionRK4`, the particle still moves down, because the grid itself moves down\n" ] } ], "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.11.6" } }, "nbformat": 4, "nbformat_minor": 2 }