{ "cells": [ { "attachments": {}, "cell_type": "markdown", "metadata": {}, "source": [ "# Delayed starts\n" ] }, { "attachments": {}, "cell_type": "markdown", "metadata": {}, "source": [ "In many applications, it is needed to 'delay' the start of particle advection. For example because particles need to be released at different times throughout an experiment. Or because particles need to be released at a constant rate from the same set of locations.\n", "\n", "This tutorial will show how this can be done. We start with importing the relevant modules.\n" ] }, { "cell_type": "code", "execution_count": 1, "metadata": {}, "outputs": [], "source": [ "from datetime import timedelta as delta\n", "\n", "import matplotlib.pyplot as plt\n", "import numpy as np\n", "import xarray as xr\n", "from IPython.display import HTML\n", "from matplotlib.animation import FuncAnimation\n", "\n", "from parcels import (\n", " AdvectionRK4,\n", " FieldSet,\n", " JITParticle,\n", " ParticleFile,\n", " ParticleSet,\n", " Variable,\n", " download_example_dataset,\n", ")" ] }, { "attachments": {}, "cell_type": "markdown", "metadata": {}, "source": [ "First import a `FieldSet` (from the Peninsula example, in this case)\n" ] }, { "cell_type": "code", "execution_count": 2, "metadata": {}, "outputs": [], "source": [ "example_dataset_folder = download_example_dataset(\"Peninsula_data\")\n", "fieldset = FieldSet.from_parcels(\n", " f\"{example_dataset_folder}/peninsula\", allow_time_extrapolation=True\n", ")" ] }, { "attachments": {}, "cell_type": "markdown", "metadata": {}, "source": [ "Now, there are two ways to delay the start of particles. Either by defining the whole `ParticleSet` at initialisation and giving each particle its own `time`. Or by using the `repeatdt` argument. We will show both options here\n" ] }, { "attachments": {}, "cell_type": "markdown", "metadata": {}, "source": [ "## Assigning each particle its own `time`\n" ] }, { "attachments": {}, "cell_type": "markdown", "metadata": {}, "source": [ "The simplest way to delaye the start of a particle is to use the `time` argument for each particle\n" ] }, { "cell_type": "code", "execution_count": 3, "metadata": {}, "outputs": [], "source": [ "npart = 10 # number of particles to be released\n", "lon = 3e3 * np.ones(npart)\n", "lat = np.linspace(3e3, 45e3, npart, dtype=np.float32)\n", "\n", "# release every particle one hour later\n", "time = np.arange(0, npart) * delta(hours=1).total_seconds()\n", "\n", "pset = ParticleSet(fieldset=fieldset, pclass=JITParticle, lon=lon, lat=lat, time=time)" ] }, { "attachments": {}, "cell_type": "markdown", "metadata": {}, "source": [ "Then we can execute the `pset` as usual\n" ] }, { "cell_type": "code", "execution_count": 4, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "INFO: Output files are stored in DelayParticle_time.zarr.\n", "100%|██████████| 86400.0/86400.0 [00:02<00:00, 35968.90it/s]\n" ] } ], "source": [ "output_file = pset.ParticleFile(name=\"DelayParticle_time.zarr\", outputdt=delta(hours=1))\n", "\n", "pset.execute(\n", " AdvectionRK4, runtime=delta(hours=24), dt=delta(minutes=5), output_file=output_file\n", ")" ] }, { "attachments": {}, "cell_type": "markdown", "metadata": {}, "source": [ "And then finally, we can show a movie of the particles. Note that the southern-most particles start to move first.\n" ] }, { "cell_type": "code", "execution_count": 5, "metadata": {}, "outputs": [], "source": [ "%%capture\n", "\n", "ds = xr.open_zarr(\"DelayParticle_time.zarr\")\n", "\n", "fig = plt.figure(figsize=(7, 5), constrained_layout=True)\n", "ax = fig.add_subplot()\n", "\n", "ax.set_ylabel(\"Meridional distance [m]\")\n", "ax.set_xlabel(\"Zonal distance [m]\")\n", "ax.set_xlim(0, 9e4)\n", "ax.set_ylim(0, 5e4)\n", "\n", "timerange = np.unique(ds[\"time\"].values[np.isfinite(ds[\"time\"])])\n", "\n", "# Indices of the data where time = 0\n", "time_id = np.where(ds[\"time\"] == timerange[0])\n", "\n", "sc = ax.scatter(ds[\"lon\"].values[time_id], ds[\"lat\"].values[time_id])\n", "\n", "t = str(timerange[0].astype(\"timedelta64[h]\"))\n", "title = ax.set_title(f\"Particles at t = {t}\")\n", "\n", "\n", "def animate(i):\n", " t = str(timerange[i].astype(\"timedelta64[h]\"))\n", " title.set_text(f\"Particles at t = {t}\")\n", "\n", " time_id = np.where(ds[\"time\"] == timerange[i])\n", " sc.set_offsets(np.c_[ds[\"lon\"].values[time_id], ds[\"lat\"].values[time_id]])\n", "\n", "\n", "anim = FuncAnimation(fig, animate, frames=len(timerange), interval=100)" ] }, { "cell_type": "code", "execution_count": 6, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\n", "\n", "\n", "\n", "\n", "\n", "