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.
This tutorial will show how this can be done. We start with importing the relevant modules.
%matplotlib inline
from parcels import FieldSet, ParticleSet, JITParticle, Variable, ParticleFile, plotTrajectoriesFile
from parcels import AdvectionRK4
import numpy as np
import xarray as xr
from datetime import timedelta as delta
import matplotlib.pyplot as plt
First import a FieldSet
(from the Peninsula example, in this case)
fieldset = FieldSet.from_parcels('Peninsula_data/peninsula', allow_time_extrapolation = True)
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
time
¶The simplest way to delaye the start of a particle is to use the time
argument for each particle
npart = 10 # number of particles to be released
lon = 3e3 * np.ones(npart)
lat = np.linspace(3e3 , 45e3, npart, dtype=np.float32)
time = np.arange(0, npart) * delta(hours=1).total_seconds() # release every particle one hour later
pset = ParticleSet(fieldset=fieldset, pclass=JITParticle, lon=lon, lat=lat, time=time)
Then we can execute the pset
as usual
output_file = pset.ParticleFile(name="DelayParticle_time.zarr", outputdt=delta(hours=1))
pset.execute(AdvectionRK4, runtime=delta(hours=24), dt=delta(minutes=5),
output_file=output_file)
INFO: Compiled ArrayJITParticleAdvectionRK4 ==> /var/folders/r2/8593q8z93kd7t4j9kbb_f7p00000gr/T/parcels-504/lib15ab640308c954e580b6f93be742baae_0.so
And then finally, we can show a movie of the particles. Note that the southern-most particles start to move first.
plotTrajectoriesFile('DelayParticle_time.zarr', mode='movie2d_notebook')
repeatdt
argument¶The second method to delay the start of particle releases is to use the repeatdt
argument when constructing a ParticleSet
. This is especially useful if you want to repeatedly release particles from the same set of locations.
npart = 10 # number of particles to be released
lon = 3e3 * np.ones(npart)
lat = np.linspace(3e3 , 45e3, npart, dtype=np.float32)
repeatdt = delta(hours=3) # release from the same set of locations every 3 hours
pset = ParticleSet(fieldset=fieldset, pclass=JITParticle, lon=lon, lat=lat, repeatdt=repeatdt)
Now we again define an output file and execute the pset
as usual.
output_file = pset.ParticleFile(name="DelayParticle_releasedt", outputdt=delta(hours=1))
pset.execute(AdvectionRK4, runtime=delta(hours=24), dt=delta(minutes=5),
output_file=output_file)
INFO: Compiled ArrayJITParticleAdvectionRK4 ==> /var/folders/r2/8593q8z93kd7t4j9kbb_f7p00000gr/T/parcels-504/libcbcec39cf73cfa7671e8798abaac733b_0.so
And we get an animation where a new particle is released every 3 hours from each start location
plotTrajectoriesFile('DelayParticle_releasedt.zarr', mode='movie2d_notebook')
Note that, if you want to if you want to at some point stop the repeatdt, the easiest implementation is to use two calls to pset.execute()
. For example, if in the above example you only want four releases of the pset, you could do the following
pset = ParticleSet(fieldset=fieldset, pclass=JITParticle, lon=lon, lat=lat, repeatdt=repeatdt)
output_file = pset.ParticleFile(name="DelayParticle_releasedt_9hrs", outputdt=delta(hours=1))
# first run for 3 * 3 hrs
pset.execute(AdvectionRK4, runtime=delta(hours=9), dt=delta(minutes=5),
output_file=output_file)
# now stop the repeated release
pset.repeatdt = None
# now continue running for the remaining 15 hours
pset.execute(AdvectionRK4, runtime=delta(hours=15), dt=delta(minutes=5),
output_file=output_file)
plotTrajectoriesFile('DelayParticle_releasedt_9hrs.zarr', mode='movie2d_notebook')
INFO: Compiled ArrayJITParticleAdvectionRK4 ==> /var/folders/r2/8593q8z93kd7t4j9kbb_f7p00000gr/T/parcels-504/libd1c17d6554e6d39a5c989beb50bcc2d3_0.so
time
in the output file¶Note that, because the outputdt
variable controls the JIT-loop, all particles are written at the same time, even when they start at a non-multiple of outputdt
.
For example, if your particles start at time=[0, 1, 2]
and outputdt=2
, then the times written (for dt=1
and endtime=4
) will be
outtime_expected = np.array([[0, 2, 4], [2, 4, np.datetime64("NaT")], [2, 4, np.datetime64("NaT")]], dtype="timedelta64[s]")
print(outtime_expected)
[[ 0 2 4] [ 2 4 'NaT'] [ 2 4 'NaT']]
outfilepath = "DelayParticle_nonmatchingtime.zarr"
pset = ParticleSet(fieldset=fieldset, pclass=JITParticle,
lat=[3e3]*3, lon=[3e3]*3, time=[0, 1, 2])
pset.execute(AdvectionRK4, endtime=4, dt=1, output_file=pset.ParticleFile(name=outfilepath, outputdt=2))
INFO: Compiled ArrayJITParticleAdvectionRK4 ==> /var/folders/r2/8593q8z93kd7t4j9kbb_f7p00000gr/T/parcels-504/lib73f1058452aa339434e56df780d0ba46_0.so
And indeed, the time
values in the NetCDF output file are as expected
outtime_infile = xr.open_dataset(outfilepath).time.values[:]
print(outtime_infile.astype('timedelta64[s]'))
assert (outtime_expected[np.isfinite(outtime_expected)] == outtime_infile[np.isfinite(outtime_infile)]).all()
[[ 0 2 4] [ 2 4 'NaT'] [ 2 4 'NaT']]
Now, for some applications, this behavior may be undesirable; for example when particles need to be analyzed at a same age (instead of at a same time). In that case, we recommend either changing outputdt
so that it is a common divisor of all start times; or doing multiple Parcels runs with subsets of the original ParticleSet
(e.g., in the example above, one run with the Particles that start at time=[0, 2]
and one with the Particle at time=[1]
). In that case, you will get two files:
for times in [[0,2], [1]]:
pset = ParticleSet(fieldset=fieldset, pclass=JITParticle,
lat=[3e3]*len(times), lon=[3e3]*len(times), time=times)
pset.execute(AdvectionRK4, endtime=4, dt=1, output_file=pset.ParticleFile(name=outfilepath, outputdt=2))
print(xr.open_dataset(outfilepath).time.values[:].astype('timedelta64[s]'))
INFO: Compiled ArrayJITParticleAdvectionRK4 ==> /var/folders/r2/8593q8z93kd7t4j9kbb_f7p00000gr/T/parcels-504/lib4a9283c7ab81279666f072b74de9d644_0.so
[[ 0 2 4] [ 2 4 'NaT']]
INFO: Compiled ArrayJITParticleAdvectionRK4 ==> /var/folders/r2/8593q8z93kd7t4j9kbb_f7p00000gr/T/parcels-504/libef5ef37165aa5febcc2e916a89504269_0.so
[[1 3 4]]
In the examples above, all particles were defined at the start of the simulation. There are use-cases, though, where it is important to be able to add particles 'on-the-fly', during the runtime of a Parcels simulation.
Unfortuantely, Parcels does not (yet) support adding new particles within a Kernel
. A workaround is to temporarily leave the execution()
call to add particles via the ParticleSet.add()
method, before continuing with execution()
.
See the example below, where we add 'mass' to a particle each timestep, based on a probablistic condition, and then split the particle once its 'mass' is larger than 5
class GrowingParticle(JITParticle):
mass = Variable('mass', initial=0)
splittime = Variable('splittime', initial=-1)
splitmass = Variable('splitmass', initial=0)
def GrowParticles(particle, fieldset, time):
if ParcelsRandom.random() < 0.25: # 25% chance per timestep for particle to grow
particle.mass += 1.
if (particle.mass >= 5.) and (particle.splittime < 0):
particle.splittime = time
particle.splitmass = particle.mass / 2.
particle.mass = particle.mass / 2.
pset = ParticleSet(fieldset=fieldset, pclass=GrowingParticle, lon=0, lat=0)
outfile = ParticleFile('growingparticles.zarr', pset, outputdt=1)
for t in range(40):
pset.execute(GrowParticles, runtime=1, dt=1, output_file=outfile)
for p in pset:
if p.splittime > 0:
pset.add(ParticleSet(fieldset=fieldset, pclass=GrowingParticle, lon=0, lat=0,
time=p.splittime, mass=p.splitmass))
p.splittime = -1 # reset splittime
INFO: Compiled ArrayGrowingParticleGrowParticles ==> /var/folders/s5/gxtkk3c12yqgd7hkt1b_s0vr0000gq/T/parcels-503/libb16839740d01629200df57019f843d0e_0.so
The 'trick' is that we place the pset.execute()
call in a for-loop, so that we leave the JIT-mode and can add Particles to the ParticleSet.
Indeed, if we plot the mass of particles as a function of time, we see that new particles are created every time a particle reaches a mass of 5.
ds = xr.open_zarr('growingparticles.zarr')
plt.plot(ds.time.values[:].astype('timedelta64[s]').T, ds.mass.T)
plt.grid()
plt.xlabel('Time')
plt.ylabel("Particle 'mass'")
plt.show()
ds.close()