User Guide

As shown in the Examples section, the typical worflow is:

  1. set up a FlowFields data structure (𝐹)
  2. set up Individuals (𝐼) with initial position 📌 and 𝐹
  3. displace 𝐼.📌 by ∫!(𝐼,𝑇) following 𝐼.𝐹 over 𝑇
  4. post-process by 𝐼.🔧 and record information in 𝐼.🔴
  5. go back to step 2 and continue if needed

The data structures for steps 1 and 2 are documented below. Steps 3 and 4 normally take place as part of ∫! which post-processes results, using 🔧, records them in 🔴, and finally updates the positions of individuals in 📌. Since 🔴 is a DataFrame, it is easily manipulated, plotted, or saved in step 4 or after the fact.

Overview

A central goal of this package is to support scientific analysis of model simulations and observations of materials and particles within the Climate System. The package scope thus includes, for example, drifting plastics in the Ocean or chemical compounds in the Atmosphere.

As a starting point, the package supports all types of gridded model output (on Arakawa C-grids) from the MIT General Circulation Model via the MeshArrays.jl package (docs found here).

By convention, IndividualDisplacements.jl expects input flow fields to be provided in a uniform fashion (see FlowFields) summarized below:

  1. normalized to grid index units (i.e. in 1/s rather than m/s units)
  2. positive towards increasing indices
  3. using the Arakawa C-grid, with u (resp v) staggered by -0.5 point in direction 1 (resp 2) from grid cell centers.

The Examples section documents various, simple methods to prepare and ingest such flow fields (time varying or not; in 2D or 3D) and derive individual displacements / trajectories from them. They cover simple grids often used in process studies, Global Ocean simulations normally done on more complex grids, plotting tools, and data formats.

For an overview of the examples, please refer to the example guide. The rest of this section is focused on the package's data structures and core functions.

Data Structures

The Individuals struct contains a FlowFields struct (incl. e.g. arrays), initial positions for the individuals, and the other elements (e.g. functions) involved in ∫!(𝐼,𝑇) as documented hereafter.

IndividualDisplacements.FlowFieldsType
abstract type FlowFields

Data structure that provide access to flow fields (gridded as arrays) which will be used to interpolate velocities to individual locations later on (once embedded in an Individuals struct).

Following the C-grid convention also used in MITgcm (https://mitgcm.readthedocs.io) flow fields are expected to be staggered as follows: grid cell i,j has its center located at i-1/2,j-1/2 while the corresponding u[i,j] (resp. `v[i,j]) is located at i-1,j-1/2 (resp. i-1/2,j-1).

Also by convention, velocity fields are expected to have been normalized to grid units (e.g. 1/s rather than m/s) before sending them to one of the supported FlowFields constructors (using either Array or MeshArray):

𝐹_Array2D(u0,u1,v0,v1,𝑇)
𝐹_Array3D(u0,u1,v0,v1,w0,w1,𝑇)
𝐹_MeshArray2D(u0,u1,v0,v1,𝑇,update_location!)
𝐹_MeshArray3D(u0,u1,v0,v1,w0,w1,𝑇,update_location!)

Using the FlowFields constructor which gets selected by the type of u0 etc. For example :

𝐹=FlowFields(u,u,v,v,0*w,1*w,[0.0,10.0])
𝐹=FlowFields(u,u,v,v,[0.0,10.0],func)

as shown in the online documentation examples.

source
IndividualDisplacements.IndividualsType
struct Individuals{T,N}
  • Data: 📌 (position), 🔴(record), 🆔 (ID), 𝑃 (FlowFields)
  • Functions: 🚄 (velocity), ∫ (integration), 🔧(postprocessing)
  • NamedTuples: 𝐷 (diagnostics), 𝑀 (metadata)

The velocity function 🚄 typically computes velocity at individual positions (📌 to start) within the specified space-time domain by interpolating gridded variables (provided via 𝑃). Individual trajectories are computed by integrating (∫) interpolated velocities through time. Normally, integration is done by calling ∫! which updates 📌 at the end and records results in 🔴 via 🔧. Ancillary data, for use in 🔧 for example, can be provided in 𝐷 and metadata stored in 𝑀.

Unicode cheatsheet:

  • 📌=\:pushpin:<tab>, 🔴=\:red_circle:<tab>, 🆔=\:id:<tab>
  • 🚄=\:bullettrain_side:<tab>, ∫=\int<tab>, 🔧=\:wrench:<tab>
  • 𝑃=\itP<tab>, 𝐷=\itD<tab>, 𝑀=\itM<tab>

Simple constructors that use FlowFields to choose adequate defaults:

  • Individuals(𝐹::𝐹_Array2D,x,y)
  • Individuals(𝐹::𝐹_Array3D,x,y,z)
  • Individuals(𝐹::𝐹_MeshArray2D,x,y,fid)
  • Individuals(𝐹::𝐹_MeshArray3D,x,y,z,fid)

Further customization is achievable via keyword constructors:

df=DataFrame( ID=[], x=[], y=[], z=[], t = [])
𝐼=Individuals{Float64,2}(📌=zeros(3,10),🆔=1:10,🔴=deepcopy(df))
𝐼=Individuals(📌=zeros(3,2),🆔=collect(1:2),🔴=deepcopy(df))

Or via the plain text (or no-unicode) constructors:

df=DataFrame( ID=[], x=[], y=[], z=[], t = [])
I=(position=zeros(3,2),ID=1:2,record=deepcopy(df))
I=Individuals(I)
source

Main Functions

∫!(𝐼,𝑇) displaces individuals 𝐼 continuously over time period 𝑇 according to velocity function 🚄, temporal integration method ∫, and post-processor 🔧 (all embedded within 𝐼).

IndividualDisplacements.∫!Function
∫!(𝐼::Individuals,𝑇::Tuple)

Displace simulated individuals continuously through space over time period 𝑇 starting from position 📌.

  • This is typically achieved by computing the cumulative integral of velocity experienced by each individual along its trajectory (∫ 🚄 dt).
  • The current default is solve(prob,Tsit5(),reltol=1e-8,abstol=1e-8) but all solver options from the OrdinaryDiffEq.jl package are available.
  • After this, ∫! is also equipped to postprocess results recorded into 🔴 via the 🔧 workflow, and the last step in ∫! consists in updating 📌 to be ready for continuing in a subsequent call to ∫!.
source
∫!(𝐼::Individuals)

Call ∫!(𝐼::Individuals,𝐼.𝑃.𝑇)

source

The velocity interpolation functions (🚄) carry out the central computation of this package – interpolating gridded flow fields to individual positions 📌. It is normally called via ∫! to integrate velocity 🚄 over a chosen time period.

  • Velocity interpolation for several array and grid types.
  • Preprocessing and postprocessing methods.
  • I/O routines to read (write) results from (to) file.

and other functionalities provided in src/compute.jl and src/data_wrangling.jl are further documented in the Tool Box section. Ingestion of trajectory data which have been collected by the Ocean Drifting Buoy Program (movie) is also supported.