LT Toolbox Tutorial: 2. Advanced Analysis with Trajectory Frames#

Welcome to the Lagrangian Trajectories Toolbox tutorial!

The LT Toolbox is a Python library dedicated to the post-processing, visualisation and analysis of Lagrangian water parcel trajectories. The toolbox offers users two structures for working with Lagrangian trajectories: Trajectory Arrays (TrajArrays) and Trajectory Frames (TrajFrames). In this tutorial, we continue exploring TrajFrames, which make use of Polars blazingly fast DataFrames to store column variables associated with trajectories (e.g. lat, lon, in-situ temperature etc.).

In this advanced tutorial, we will learn how to:

  • Add new attribute variables using the Polars expressions syntax to your TrajFrame.

  • Compute binned statistics, such as the Lagrangian probability, and grouped statistics, such as meridional heat transport, from Lagrangian trajectories.

  • Map probabilities and properties determined from Lagrangian trajectories.

To learn about the basics of working with TrajFrames, including simple filtering and using datetimes, users should see LT Toolbox Tutorial: 1. Getting Started with Trajectory Frames.

Getting Started#

Let us begin by importing the relevant packages we’ll need to get started with the LT Toolbox.

Note: Since lt_toolbox is still undergoing unit testing, the package is not yet available on PyPi, we use pip to install a local development version.

[1]:
# Import relevant Python packages:
import xarray as xr
import numpy as np
import polars as pl
import matplotlib.pyplot as plt

# Following pip installation as shown on the LT Toolbox github:
import lt_toolbox as ltt

Storing Trajectory Data#

To explore the functionality of the LT Toolbox, we will use example output from a Lagrangian particle tracking experiment using the eddy-rich ORCA0083-GO8p7 JRA55 ocean sea-ice hindcast (Megann et al. 2022). Trajectories were advected forwards-in-time from the full-depth northward inflows across the Overturning in the Subpolar North Atlantic (OSNAP) East array.

We have provided some example trajectory output for the January 1995 initialisation in the .parquet file format.

Below we load the output .parquet file as a DataFrame with Polars, before creating a TrajFrame, traj.

[2]:
# Defining filepath to our example ORCA0083-GO8p7 output trajectory file:
traj_filepath = "./data/ORCA0083-GO8p7_JRA55_SPNA_1995_example.parquet"

# Open output .parquet file as a DataFrame.
dataset = pl.read_parquet(traj_filepath, use_pyarrow=True)

# Create a TrajFrame from the DataFrame:
traj = ltt.TrajFrame(source=dataset, condense=True)

# Here we see how condensing a TrajFrame will store all observations in
# lists per trajectory (rows) instead of in a row per observation.
traj
[2]:
<TrajFrame object>

----- Trajectory DataFrame -----
Trajectories: 2500
Variables: ['id', 'x', 'y', 'z', 'subvol', 'time', 'boxface', 'thetao', 'so', 'mask']
<bound method DataFrame.glimpse of shape: (2_500, 10)
┌───────┬────────────┬────────────┬────────────┬───┬───────────┬───────────┬───────────┬───────────┐
│ id    ┆ x          ┆ y          ┆ z          ┆ … ┆ boxface   ┆ thetao    ┆ so        ┆ mask      │
│ ---   ┆ ---        ┆ ---        ┆ ---        ┆   ┆ ---       ┆ ---       ┆ ---       ┆ ---       │
│ i64   ┆ list[f64]  ┆ list[f64]  ┆ list[f64]  ┆   ┆ list[i64] ┆ list[f64] ┆ list[f64] ┆ list[f64] │
╞═══════╪════════════╪════════════╪════════════╪═══╪═══════════╪═══════════╪═══════════╪═══════════╡
│ 1     ┆ [3085.0,   ┆ [2888.5,   ┆ [46.5,     ┆ … ┆ [1, 4, …  ┆ [4.35,    ┆ [35.22,   ┆ [1.0,     │
│       ┆ 3085.0, …  ┆ 2888.0, …  ┆ 46.42, …   ┆   ┆ 4]        ┆ 4.34, …   ┆ 35.22, …  ┆ 1.0, …    │
│       ┆ 2840.83]   ┆ 2793.0]    ┆ 48.03]     ┆   ┆           ┆ 3.81]     ┆ 35.16]    ┆ 0.5]      │
│ 12    ┆ [3085.0,   ┆ [2888.83,  ┆ [48.83,    ┆ … ┆ [1, 4, …  ┆ [4.02,    ┆ [35.25,   ┆ [1.0,     │
│       ┆ 3085.34, … ┆ 2888.0, …  ┆ 48.14, …   ┆   ┆ 4]        ┆ 4.01, …   ┆ 35.25, …  ┆ 1.0, …    │
│       ┆ 2847.51]   ┆ 2797.0]    ┆ 58.28]     ┆   ┆           ┆ 2.84]     ┆ 35.16]    ┆ 0.5]      │
│ 15    ┆ [3085.0,   ┆ [2888.12,  ┆ [49.87,    ┆ … ┆ [1, 6, …  ┆ [4.0,     ┆ [35.25,   ┆ [1.0,     │
│       ┆ 3085.98, … ┆ 2888.07, … ┆ 49.0, …    ┆   ┆ 1]        ┆ 4.09, …   ┆ 35.25, …  ┆ 1.0, …    │
│       ┆ 2844.0]    ┆ 2794.02]   ┆ 46.57]     ┆   ┆           ┆ 3.56]     ┆ 35.17]    ┆ 0.0]      │
│ 21    ┆ [3085.0,   ┆ [2888.37,  ┆ [49.37,    ┆ … ┆ [1, 6, …  ┆ [4.0,     ┆ [35.25,   ┆ [1.0,     │
│       ┆ 3085.36, … ┆ 2888.3, …  ┆ 49.0, …    ┆   ┆ 4]        ┆ 4.1, …    ┆ 35.25, …  ┆ 1.0, …    │
│       ┆ 2842.74]   ┆ 2794.0]    ┆ 52.38]     ┆   ┆           ┆ 3.35]     ┆ 35.18]    ┆ 0.5]      │
│ 26    ┆ [3085.0,   ┆ [2888.62,  ┆ [49.12,    ┆ … ┆ [1, 6, …  ┆ [4.0,     ┆ [35.25,   ┆ [1.0,     │
│       ┆ 3085.11, … ┆ 2888.58, … ┆ 49.0, …    ┆   ┆ 4]        ┆ 4.1, …    ┆ 35.25, …  ┆ 1.0, …    │
│       ┆ 2843.97]   ┆ 2794.0]    ┆ 50.7]      ┆   ┆           ┆ 3.52]     ┆ 35.19]    ┆ 0.5]      │
│ …     ┆ …          ┆ …          ┆ …          ┆ … ┆ …         ┆ …         ┆ …         ┆ …         │
│ 28473 ┆ [3311.75,  ┆ [2865.0,   ┆ [37.75,    ┆ … ┆ [3, 5, …  ┆ [10.22,   ┆ [35.54,   ┆ [1.0,     │
│       ┆ 3311.52, … ┆ 2865.4, …  ┆ 38.0, …    ┆   ┆ 4]        ┆ 10.2, …   ┆ 35.54, …  ┆ 1.0, …    │
│       ┆ 2844.18]   ┆ 2795.0]    ┆ 54.49]     ┆   ┆           ┆ 3.21]     ┆ 35.17]    ┆ 0.5]      │
│ 28512 ┆ [3312.83,  ┆ [2864.0,   ┆ [24.5,     ┆ … ┆ [3, 6, …  ┆ [10.27,   ┆ [35.53,   ┆ [1.0,     │
│       ┆ 3312.65, … ┆ 2864.56, … ┆ 24.0, …    ┆   ┆ 1]        ┆ 9.97, …   ┆ 35.53, …  ┆ 1.0, …    │
│       ┆ 2840.0]    ┆ 2792.97]   ┆ 40.77]     ┆   ┆           ┆ 4.67]     ┆ 35.11]    ┆ 0.0]      │
│ 28549 ┆ [3312.5,   ┆ [2864.0,   ┆ [35.5,     ┆ … ┆ [3, 2, …  ┆ [10.26,   ┆ [35.53,   ┆ [1.0,     │
│       ┆ 3312.0, …  ┆ 2864.52, … ┆ 35.47, …   ┆   ┆ 4]        ┆ 10.23, …  ┆ 35.53, …  ┆ 1.0, …    │
│       ┆ 2844.09]   ┆ 2795.0]    ┆ 47.4]      ┆   ┆           ┆ 3.57]     ┆ 35.18]    ┆ 0.5]      │
│ 28956 ┆ [3322.5,   ┆ [2862.0,   ┆ [11.5,     ┆ … ┆ [3, 1, …  ┆ [9.61,    ┆ [35.35,   ┆ [1.0,     │
│       ┆ 3323.0, …  ┆ 2862.47, … ┆ 11.65, …   ┆   ┆ 4]        ┆ 8.72, …   ┆ 35.33, …  ┆ 1.0, …    │
│       ┆ 2842.39]   ┆ 2794.0]    ┆ 48.38]     ┆   ┆           ┆ 3.69]     ┆ 35.18]    ┆ 0.5]      │
│ 28990 ┆ [3323.5,   ┆ [2862.0,   ┆ [18.5,     ┆ … ┆ [3, 1, …  ┆ [9.54,    ┆ [35.32,   ┆ [1.0,     │
│       ┆ 3324.0, …  ┆ 2862.27, … ┆ 18.56, …   ┆   ┆ 4]        ┆ 8.65, …   ┆ 35.32, …  ┆ 1.0, …    │
│       ┆ 2845.52]   ┆ 2795.0]    ┆ 49.83]     ┆   ┆           ┆ 3.54]     ┆ 35.19]    ┆ 0.5]      │
└───────┴────────────┴────────────┴────────────┴───┴───────────┴───────────┴───────────┴───────────┘>

Transforming x,y,z variables to lon,lat,depth in a TrajFrame.#

In our example TrajFrame, we have stored the positions of each trajectory as indices referring to the original numerical model grid (eORCA12). However, we can easily transform these model coordinates to geographical coordinates using the .transform_trajectory_coords() method which uses bi-linear interpolation to determine each trajectories geographical position.

[3]:
# Importing eORCA12 ocean model grid from a NetCDF file using xarray:
ds_grid = xr.open_dataset('/home/snapdragon/HadGEM3-GC31-MM/Proj_Future_Pathways/src/Software/lt_toolbox/docs/tutorials/data/ORCA0083-GO8p7_JRA55_model_grid.nc')

ds_grid
[3]:
<xarray.Dataset> Size: 125MB
Dimensions:  (y: 3606, x: 4322, z: 75)
Dimensions without coordinates: y, x, z
Data variables:
    nav_lat  (y, x) float32 62MB ...
    nav_lon  (y, x) float32 62MB ...
    nav_lev  (z) float32 300B ...
Attributes:
    file_name:  mesh_mask.nc
    TimeStamp:  30/09/2016 08:33:26 +0000
[4]:
# Transform the TrajFrame to geographical coordinates and store as new TrajFrame object:
traj_geo = traj.transform_trajectory_coords(lon=ds_grid.nav_lon.values,
                                            lat=ds_grid.nav_lat.values,
                                            depth=ds_grid.nav_lev.values,
                                            )

traj_geo
[4]:
<TrajFrame object>

----- Trajectory DataFrame -----
Trajectories: 2500
Variables: ['id', 'lon', 'lat', 'depth', 'subvol', 'time', 'boxface', 'thetao', 'so', 'mask']
<bound method DataFrame.glimpse of shape: (2_500, 10)
┌───────┬────────────┬────────────┬────────────┬───┬───────────┬───────────┬───────────┬───────────┐
│ id    ┆ lon        ┆ lat        ┆ depth      ┆ … ┆ boxface   ┆ thetao    ┆ so        ┆ mask      │
│ ---   ┆ ---        ┆ ---        ┆ ---        ┆   ┆ ---       ┆ ---       ┆ ---       ┆ ---       │
│ i64   ┆ list[f64]  ┆ list[f64]  ┆ list[f64]  ┆   ┆ list[i64] ┆ list[f64] ┆ list[f64] ┆ list[f64] │
╞═══════╪════════════╪════════════╪════════════╪═══╪═══════════╪═══════════╪═══════════╪═══════════╡
│ 1     ┆ [-30.38194 ┆ [58.799782 ┆ [1098.9227 ┆ … ┆ [1, 4, …  ┆ [4.35,    ┆ [35.22,   ┆ [1.0,     │
│       ┆ 3, -30.380 ┆ ,          ┆ 29, 1090.4 ┆   ┆ 4]        ┆ 4.34, …   ┆ 35.22, …  ┆ 1.0, …    │
│       ┆ 463, …     ┆ 58.774166, ┆ 31772, …   ┆   ┆           ┆ 3.81]     ┆ 35.16]    ┆ 0.5]      │
│       ┆ -51.5…     ┆ … 52.9276… ┆ 126…       ┆   ┆           ┆           ┆           ┆           │
│ 12    ┆ [-30.38292 ┆ [58.816688 ┆ [1366.7193 ┆ … ┆ [1, 4, …  ┆ [4.02,    ┆ [35.25,   ┆ [1.0,     │
│       ┆ , -30.3485 ┆ ,          ┆ 18, 1282.8 ┆   ┆ 4]        ┆ 4.01, …   ┆ 35.25, …  ┆ 1.0, …    │
│       ┆ 98, …      ┆ 58.774656, ┆ 73621, …   ┆   ┆           ┆ 2.84]     ┆ 35.16]    ┆ 0.5]      │
│       ┆ -51.02…    ┆ … 53.1688… ┆ 282…       ┆   ┆           ┆           ┆           ┆           │
│ 15    ┆ [-30.38081 ┆ [58.780314 ┆ [1499.5953 ┆ … ┆ [1, 6, …  ┆ [4.0,     ┆ [35.25,   ┆ [1.0,     │
│       ┆ 8, -30.288 ┆ ,          ┆ 77, 1387.3 ┆   ┆ 1]        ┆ 4.09, …   ┆ 35.25, …  ┆ 1.0, …    │
│       ┆ 821, …     ┆ 58.779166, ┆ 76953, …   ┆   ┆           ┆ 3.56]     ┆ 35.17]    ┆ 0.0]      │
│       ┆ -51.3…     ┆ … 52.9971… ┆ 110…       ┆   ┆           ┆           ┆           ┆           │
│ 21    ┆ [-30.38155 ┆ [58.793122 ┆ [1435.1020 ┆ … ┆ [1, 6, …  ┆ [4.0,     ┆ [35.25,   ┆ [1.0,     │
│       ┆ 8, -30.347 ┆ ,          ┆ 3, 1387.37 ┆   ┆ 4]        ┆ 4.1, …    ┆ 35.25, …  ┆ 1.0, …    │
│       ┆ 605, …     ┆ 58.790056, ┆ 6953, …    ┆   ┆           ┆ 3.35]     ┆ 35.18]    ┆ 0.5]      │
│       ┆ -51.4…     ┆ … 52.9893… ┆ 1852…      ┆   ┆           ┆           ┆           ┆           │
│ 26    ┆ [-30.38229 ┆ [58.80593, ┆ [1402.8553 ┆ … ┆ [1, 6, …  ┆ [4.0,     ┆ [35.25,   ┆ [1.0,     │
│       ┆ 8, -30.371 ┆ 58.80404,  ┆ 56, 1387.3 ┆   ┆ 4]        ┆ 4.1, …    ┆ 35.25, …  ┆ 1.0, …    │
│       ┆ 867, …     ┆ …          ┆ 76953, …   ┆   ┆           ┆ 3.52]     ┆ 35.19]    ┆ 0.5]      │
│       ┆ -51.3…     ┆ 52.995958… ┆ 161…       ┆   ┆           ┆           ┆           ┆           │
│ …     ┆ …          ┆ …          ┆ …          ┆ … ┆ …         ┆ …         ┆ …         ┆ …         │
│ 28473 ┆ [-9.324242 ┆ [57.098342 ┆ [446.16767 ┆ … ┆ [3, 5, …  ┆ [10.22,   ┆ [35.54,   ┆ [1.0,     │
│       ┆ ,          ┆ ,          ┆ 1,         ┆   ┆ 4]        ┆ 10.2, …   ┆ 35.54, …  ┆ 1.0, …    │
│       ┆ -9.341056, ┆ 57.120076, ┆ 457.62561, ┆   ┆           ┆ 3.21]     ┆ 35.17]    ┆ 0.5]      │
│       ┆ … -51.293… ┆ … 53.0483… ┆ … 2180.1…  ┆   ┆           ┆           ┆           ┆           │
│ 28512 ┆ [-9.235759 ┆ [57.041115 ┆ [114.01514 ┆ … ┆ [3, 6, …  ┆ [10.27,   ┆ [35.53,   ┆ [1.0,     │
│       ┆ ,          ┆ ,          ┆ 1, 108.030 ┆   ┆ 1]        ┆ 9.97, …   ┆ 35.53, …  ┆ 1.0, …    │
│       ┆ -9.246296, ┆ 57.070754, ┆ 281, …     ┆   ┆           ┆ 4.67]     ┆ 35.11]    ┆ 0.0]      │
│       ┆ … -51.634… ┆ … 52.9216… ┆ 613.5…     ┆   ┆           ┆           ┆           ┆           │
│ 28549 ┆ [-9.266052 ┆ [57.042982 ┆ [352.27565 ┆ … ┆ [3, 2, …  ┆ [10.26,   ┆ [35.53,   ┆ [1.0,     │
│       ┆ ,          ┆ ,          ┆ ,          ┆   ┆ 4]        ┆ 10.23, …  ┆ 35.53, …  ┆ 1.0, …    │
│       ┆ -9.306412, ┆ 57.072393, ┆ 351.17088, ┆   ┆           ┆ 3.57]     ┆ 35.18]    ┆ 0.5]      │
│       ┆ … -51.301… ┆ … 53.0478… ┆ … 1197.53… ┆   ┆           ┆           ┆           ┆           │
│ 28956 ┆ [-8.371674 ┆ [56.882877 ┆ [21.093709 ┆ … ┆ [3, 1, …  ┆ [9.61,    ┆ [35.35,   ┆ [1.0,     │
│       ┆ ,          ┆ ,          ┆ ,          ┆   ┆ 4]        ┆ 8.72, …   ┆ 35.33, …  ┆ 1.0, …    │
│       ┆ -8.320692, ┆ 56.903886, ┆ 21.592881, ┆   ┆           ┆ 3.69]     ┆ 35.18]    ┆ 0.5]      │
│       ┆ … -51.438… ┆ … 52.9874… ┆ … 1312.03… ┆   ┆           ┆           ┆           ┆           │
│ 28990 ┆ [-8.280212 ┆ [56.876928 ┆ [57.481737 ┆ … ┆ [3, 1, …  ┆ [9.54,    ┆ [35.32,   ┆ [1.0,     │
│       ┆ ,          ┆ ,          ┆ ,          ┆   ┆ 4]        ┆ 8.65, …   ┆ 35.32, …  ┆ 1.0, …    │
│       ┆ -8.231481, ┆ 56.887714, ┆ 57.917469, ┆   ┆           ┆ 3.54]     ┆ 35.19]    ┆ 0.5]      │
│       ┆ … -51.178… ┆ … 53.0555… ┆ … 1494.43… ┆   ┆           ┆           ┆           ┆           │
└───────┴────────────┴────────────┴────────────┴───┴───────────┴───────────┴───────────┴───────────┘>

Adding Attribute Variables with List Expressions#

Often when we want to add a new attribute variables to our TrajFrame, we would like to compute a complex function using existing attribute variables, for example potential density from temperature and salinity recorded along-stream. Below we show how we can do this using Polars expressions syntax and list operations.

[5]:
# Import the density module from Gibbs Seawater Toolbox (gsw):
import gsw.density as density

# Compute potential density anomaly (sigma0) from conservative temperature
# and absolute salinity recorded along our trajectories:
traj_geo = traj_geo.add_variable(name='sigma0',
                                 expr=density.sigma0(CT=pl.col('thetao'), SA=pl.col('so')),
                                 list_expr=True,
                                 )

traj_geo
[5]:
<TrajFrame object>

----- Trajectory DataFrame -----
Trajectories: 2500
Variables: ['id', 'lon', 'lat', 'depth', 'subvol', 'time', 'boxface', 'thetao', 'so', 'mask', 'sigma0']
<bound method DataFrame.glimpse of shape: (2_500, 11)
┌───────┬────────────┬────────────┬────────────┬───┬───────────┬───────────┬───────────┬───────────┐
│ id    ┆ lon        ┆ lat        ┆ depth      ┆ … ┆ thetao    ┆ so        ┆ mask      ┆ sigma0    │
│ ---   ┆ ---        ┆ ---        ┆ ---        ┆   ┆ ---       ┆ ---       ┆ ---       ┆ ---       │
│ i64   ┆ list[f64]  ┆ list[f64]  ┆ list[f64]  ┆   ┆ list[f64] ┆ list[f64] ┆ list[f64] ┆ list[f64] │
╞═══════╪════════════╪════════════╪════════════╪═══╪═══════════╪═══════════╪═══════════╪═══════════╡
│ 1     ┆ [-30.38194 ┆ [58.799782 ┆ [1098.9227 ┆ … ┆ [4.35,    ┆ [35.22,   ┆ [1.0,     ┆ [27.79184 │
│       ┆ 3, -30.380 ┆ ,          ┆ 29, 1090.4 ┆   ┆ 4.34, …   ┆ 35.22, …  ┆ 1.0, …    ┆ 6, 27.792 │
│       ┆ 463, …     ┆ 58.774166, ┆ 31772, …   ┆   ┆ 3.81]     ┆ 35.16]    ┆ 0.5]      ┆ 939, …    │
│       ┆ -51.5…     ┆ … 52.9276… ┆ 126…       ┆   ┆           ┆           ┆           ┆ 27.8017…  │
│ 12    ┆ [-30.38292 ┆ [58.816688 ┆ [1366.7193 ┆ … ┆ [4.02,    ┆ [35.25,   ┆ [1.0,     ┆ [27.851,  │
│       ┆ , -30.3485 ┆ ,          ┆ 18, 1282.8 ┆   ┆ 4.01, …   ┆ 35.25, …  ┆ 1.0, …    ┆ 27.852055 │
│       ┆ 98, …      ┆ 58.774656, ┆ 73621, …   ┆   ┆ 2.84]     ┆ 35.16]    ┆ 0.5]      ┆ , … 27.89 │
│       ┆ -51.02…    ┆ … 53.1688… ┆ 282…       ┆   ┆           ┆           ┆           ┆ 5688]     │
│ 15    ┆ [-30.38081 ┆ [58.780314 ┆ [1499.5953 ┆ … ┆ [4.0,     ┆ [35.25,   ┆ [1.0,     ┆ [27.85310 │
│       ┆ 8, -30.288 ┆ ,          ┆ 77, 1387.3 ┆   ┆ 4.09, …   ┆ 35.25, …  ┆ 1.0, …    ┆ 8, 27.843 │
│       ┆ 821, …     ┆ 58.779166, ┆ 76953, …   ┆   ┆ 3.56]     ┆ 35.17]    ┆ 0.0]      ┆ 584, …    │
│       ┆ -51.3…     ┆ … 52.9971… ┆ 110…       ┆   ┆           ┆           ┆           ┆ 27.8349…  │
│ 21    ┆ [-30.38155 ┆ [58.793122 ┆ [1435.1020 ┆ … ┆ [4.0,     ┆ [35.25,   ┆ [1.0,     ┆ [27.85310 │
│       ┆ 8, -30.347 ┆ ,          ┆ 3, 1387.37 ┆   ┆ 4.1, …    ┆ 35.25, …  ┆ 1.0, …    ┆ 8,        │
│       ┆ 605, …     ┆ 58.790056, ┆ 6953, …    ┆   ┆ 3.35]     ┆ 35.18]    ┆ 0.5]      ┆ 27.84252, │
│       ┆ -51.4…     ┆ … 52.9893… ┆ 1852…      ┆   ┆           ┆           ┆           ┆ …         │
│       ┆            ┆            ┆            ┆   ┆           ┆           ┆           ┆ 27.86354… │
│ 26    ┆ [-30.38229 ┆ [58.80593, ┆ [1402.8553 ┆ … ┆ [4.0,     ┆ [35.25,   ┆ [1.0,     ┆ [27.85310 │
│       ┆ 8, -30.371 ┆ 58.80404,  ┆ 56, 1387.3 ┆   ┆ 4.1, …    ┆ 35.25, …  ┆ 1.0, …    ┆ 8,        │
│       ┆ 867, …     ┆ …          ┆ 76953, …   ┆   ┆ 3.52]     ┆ 35.19]    ┆ 0.5]      ┆ 27.84252, │
│       ┆ -51.3…     ┆ 52.995958… ┆ 161…       ┆   ┆           ┆           ┆           ┆ …         │
│       ┆            ┆            ┆            ┆   ┆           ┆           ┆           ┆ 27.85475… │
│ …     ┆ …          ┆ …          ┆ …          ┆ … ┆ …         ┆ …         ┆ …         ┆ …         │
│ 28473 ┆ [-9.324242 ┆ [57.098342 ┆ [446.16767 ┆ … ┆ [10.22,   ┆ [35.54,   ┆ [1.0,     ┆ [27.20444 │
│       ┆ ,          ┆ ,          ┆ 1,         ┆   ┆ 10.2, …   ┆ 35.54, …  ┆ 1.0, …    ┆ 9, 27.207 │
│       ┆ -9.341056, ┆ 57.120076, ┆ 457.62561, ┆   ┆ 3.21]     ┆ 35.17]    ┆ 0.5]      ┆ 934, …    │
│       ┆ … -51.293… ┆ … 53.0483… ┆ … 2180.1…  ┆   ┆           ┆           ┆           ┆ 27.8691…  │
│ 28512 ┆ [-9.235759 ┆ [57.041115 ┆ [114.01514 ┆ … ┆ [10.27,   ┆ [35.53,   ┆ [1.0,     ┆ [27.18798 │
│       ┆ ,          ┆ ,          ┆ 1, 108.030 ┆   ┆ 9.97, …   ┆ 35.53, …  ┆ 1.0, …    ┆ 6,        │
│       ┆ -9.246296, ┆ 57.070754, ┆ 281, …     ┆   ┆ 4.67]     ┆ 35.11]    ┆ 0.0]      ┆ 27.23998, │
│       ┆ … -51.634… ┆ … 52.9216… ┆ 613.5…     ┆   ┆           ┆           ┆           ┆ …         │
│       ┆            ┆            ┆            ┆   ┆           ┆           ┆           ┆ 27.66951… │
│ 28549 ┆ [-9.266052 ┆ [57.042982 ┆ [352.27565 ┆ … ┆ [10.26,   ┆ [35.53,   ┆ [1.0,     ┆ [27.18973 │
│       ┆ ,          ┆ ,          ┆ ,          ┆   ┆ 10.23, …  ┆ 35.53, …  ┆ 1.0, …    ┆ 3, 27.194 │
│       ┆ -9.306412, ┆ 57.072393, ┆ 351.17088, ┆   ┆ 3.57]     ┆ 35.18]    ┆ 0.5]      ┆ 971, …    │
│       ┆ … -51.301… ┆ … 53.0478… ┆ … 1197.53… ┆   ┆           ┆           ┆           ┆ 27.8418…  │
│ 28956 ┆ [-8.371674 ┆ [56.882877 ┆ [21.093709 ┆ … ┆ [9.61,    ┆ [35.35,   ┆ [1.0,     ┆ [27.16166 │
│       ┆ ,          ┆ ,          ┆ ,          ┆   ┆ 8.72, …   ┆ 35.33, …  ┆ 1.0, …    ┆ 7, 27.291 │
│       ┆ -8.320692, ┆ 56.903886, ┆ 21.592881, ┆   ┆ 3.69]     ┆ 35.18]    ┆ 0.5]      ┆ 193, …    │
│       ┆ … -51.438… ┆ … 52.9874… ┆ … 1312.03… ┆   ┆           ┆           ┆           ┆ 27.8297…  │
│ 28990 ┆ [-8.280212 ┆ [56.876928 ┆ [57.481737 ┆ … ┆ [9.54,    ┆ [35.32,   ┆ [1.0,     ┆ [27.15012 │
│       ┆ ,          ┆ ,          ┆ ,          ┆   ┆ 8.65, …   ┆ 35.32, …  ┆ 1.0, …    ┆ 1, 27.294 │
│       ┆ -8.231481, ┆ 56.887714, ┆ 57.917469, ┆   ┆ 3.54]     ┆ 35.19]    ┆ 0.5]      ┆ 473, …    │
│       ┆ … -51.178… ┆ … 53.0555… ┆ … 1494.43… ┆   ┆           ┆           ┆           ┆ 27.8527…  │
└───────┴────────────┴────────────┴────────────┴───┴───────────┴───────────┴───────────┴───────────┘>

Computing Grouped Expressions#

In addition to the more basic computations demonstrated in Tutorial 2, the LT-Toolbox also includes several more advanced computation methods which return n-dimensional arrays as outputs since they do not conform to the tabular structure of our Lagrangian trajectory data. To help manage the resulting n-dimensional arrays produced, a TrajFrame also stores these in an xarray DataSet accessible via the summary_data attribute.

Below we show a simple example of a grouped expression where we group over the water parcel ID variable in an uncondensed TrajFrame. Note, that this example is for illustrative purposes only since we could store this output in our condensed TrajFrame.

A more common example of grouped expressions is computing an aggregated statistic over groups of water parcel trajectories; for example, we could compute the average transit time grouped over start time if we had released parcels every month for a year.

[6]:
# Create an uncondensed TrajFrame from our example DataFrame:
traj = ltt.TrajFrame(source=dataset)

# Define physical constants:
rho0 = 1025  # Reference seawater density in kg/m^3
Cp = 3992  # Seawater specific heat capacity in J/kg/K

# Define Polars expression to compute the change in heat
# transport of each water parcel trajectory:
expr = rho0*Cp*pl.col('subvol').first()*(pl.col('thetao').last() - pl.col('thetao').first())

# Computing the change in heat transport of each water parcel trajectory:
traj = traj.compute_grouped_expr(group='id', expr=expr, alias='dMHT', append=False)

# Display the summary data of the TrajFrame:
traj.summary_data
[6]:
<xarray.Dataset> Size: 40kB
Dimensions:  (id: 2500)
Coordinates:
  * id       (id) int64 20kB 1 12 15 21 26 32 ... 28473 28512 28549 28956 28990
Data variables:
    dMHT     (id) float64 20kB -5.24e+08 -5.837e+09 ... -1.32e+10 -1.097e+10

Computing 1-dimensional Binned Statistics#

A further example of an advanced computation is to compute 1-dimensional binned statistics or histograms using the positions or properties sampled along water parcel trajectories.

Below we demonstrate how to compute the sum of water parcel volume transport binned according to their initial longitude at the point of release along the Overturning in the Subpolar North Atlantic Program (OSNAP) array in the subpolar North Atlantic.

We will visualise our results with xarray’s in-built plotting wrapper of the matplotlib library.

[7]:
# Define discrete longitude bins:
bin_lon = np.arange(-33, 0, 1).tolist()

# Compute 1-D summed volume transport in discrete longitude bins:
traj_geo = (traj_geo
            # Add volume transport [Sv] per parcel as an attribute variable:
            .add_variable(name='subvol_ini',
                          expr=(pl.col('subvol').list.first() / 1E6),
                          list_expr=False,
                          )
            # Add the initial longitude of each parcel as an attribute variable:
            .add_variable(name='lon_ini',
                          expr=pl.col('lon').list.first(),
                          list_expr=False,
                          )
            # Compute volume transport in discrete longitude bins:
            .compute_binned_statistic_1d(var='lon_ini',
                                         values='subvol_ini',
                                         statistic='sum',
                                         bin_breaks=bin_lon,
                                         alias='subvol_x_ini',
                                         group=None, # No group-by required.
                                         append=True, # Allow append to existing summary DataSet.
                                         )
            )

# Display the summary data of the TrajFrame:
traj_geo.summary_data
[7]:
<xarray.Dataset> Size: 368B
Dimensions:       (lon_ini: 23)
Coordinates:
  * lon_ini       (lon_ini) float64 184B -30.5 -29.5 -28.5 ... -10.5 -9.5 -8.5
Data variables:
    subvol_x_ini  (lon_ini) float64 184B 0.01537 0.8793 ... 0.02333 0.0009916
[8]:
# Plotting the summed volume transport in discrete longitude bins:
traj_geo.summary_data.subvol_x_ini.plot()
[8]:
[<matplotlib.lines.Line2D at 0x7f2850afb4c0>]
../_images/tutorials_tutorial_2_advanced_analysis_with_TrajFrames_14_1.png

Computing 2-dimensional Binned Statistics#

We next extend our computation to a 2-dimensional binned statistic using two different positions or properties sampled along water parcel trajectories.

Below we demonstrate how to compute the average potential density recorded along-stream binned according to the longitude and latitude positions of our trajectories in the subpolar North Atlantic.

We will again visualise our results with xarray’s in-built plotting wrapper of the matplotlib library.

[9]:
# Define discrete depth bins:
bin_lon = np.arange(-70, 0, 0.25).tolist()
bin_lat = np.arange(50, 70, 0.25).tolist()

# Compute 2-D mean potential density in discrete longitude-latitude bins:
traj_geo = (traj_geo
            # Compute volume transport in discrete longitude-latitude bins:
            .compute_binned_statistic_2d(var_x='lat',
                                         var_y='lon',
                                         values='sigma0',
                                         statistic='mean',
                                         bin_breaks=[bin_lat, bin_lon],
                                         alias='sigma0_xy',
                                         group=None, # No group-by required!
                                         append=True, # Allow append to existing summary DataSet.
                                         )
            )

# Display the summary data of the TrajFrame:
traj_geo.summary_data
[9]:
<xarray.Dataset> Size: 140kB
Dimensions:       (lon_ini: 23, lat: 69, lon: 248)
Coordinates:
  * lon_ini       (lon_ini) float64 184B -30.5 -29.5 -28.5 ... -10.5 -9.5 -8.5
  * lat           (lat) float64 552B 51.88 52.12 52.38 ... 68.38 68.62 68.88
  * lon           (lon) float64 2kB -66.88 -66.62 -66.38 ... -5.375 -5.125
Data variables:
    subvol_x_ini  (lon_ini) float64 184B 0.01537 0.8793 ... 0.02333 0.0009916
    sigma0_xy     (lat, lon) float64 137kB nan nan nan nan ... nan nan nan nan
[10]:
# Plotting the mean potential density along-stream in discrete 2-D longitude-latitude bins:
traj_geo.summary_data.sigma0_xy.plot()
[10]:
<matplotlib.collections.QuadMesh at 0x7f28500f66e0>
../_images/tutorials_tutorial_2_advanced_analysis_with_TrajFrames_17_1.png

Computing Lagrangian Probability Statistics#

The most common 2-dimensional binned statistic in Lagrangian oceanography is the probability density map.

Below we demonstrate how to compute two types of Lagrangian probability density map for our exmaple trajectories circulating around the North Atlantic Subpolar Gyre. The two types of Lagrangian probability are outlined in the excellent review paper by van Sebille et al. (2018):

  • ‘pos’ - What is the probability that any given observation of all water parcel trajectories will be found in a discrete longitude-latitude bin:

\(p_{pos} = \frac{N_{pos}(j,i)}{\sum_j \sum_i N_{pos}}\)

  • ‘traj’ - What is the probability that any given water parcel will enter a discrete longitude-latitude bin at least once during its trajectory:

\(p_{traj} = \frac{N_{traj}(j,i)}{\sum_j \sum_i N_{traj}}\)

We will visualise our results with xarray’s in-built plotting wrapper of the matplotlib library.

[11]:
# Compute Lagrangian probability in discrete longitude-latitude bins:
traj_geo.summary_data = xr.Dataset()
traj_geo = (traj_geo
            # Compute volume transport in discrete longitude-latitude bins:
            .compute_probability(bin_res=0.25, # Bin resolution in degrees
                                 prob_type='pos', # Type of Lagrangian probability density.
                                 group=None, # No group-by required.
                                 append=False, # Replace existing summary DataSet.
                                )
            )

# Display the summary data of the TrajFrame:
traj_geo.summary_data
[11]:
<xarray.Dataset> Size: 139kB
Dimensions:      (lat: 69, lon: 248)
Coordinates:
  * lat          (lat) float64 552B 51.88 52.12 52.38 ... 68.38 68.62 68.88
  * lon          (lon) float64 2kB -66.88 -66.62 -66.38 ... -5.625 -5.375 -5.125
Data variables:
    probability  (lat, lon) float64 137kB nan nan nan nan ... nan nan nan nan
[12]:
# Plotting the Lagrangian probability in discrete 2-D longitude-latitude bins:
traj_geo.summary_data.probability.plot()
[12]:
<matplotlib.collections.QuadMesh at 0x7f283d017430>
../_images/tutorials_tutorial_2_advanced_analysis_with_TrajFrames_20_1.png
[13]:
# Compute Lagrangian probability in discrete longitude-latitude bins:
traj_geo.summary_data = xr.Dataset()
traj_geo = (traj_geo
            # Compute volume transport in discrete longitude-latitude bins:
            .compute_probability(bin_res=0.25, # Bin resolution in degrees
                                 prob_type='traj', # Type of Lagrangian probability density
                                 group=None, # No group-by required.
                                 append=False, # Replace existing summary DataSet.
                                )
            )

# Display the summary data of the TrajFrame:
traj_geo.summary_data
[13]:
<xarray.Dataset> Size: 139kB
Dimensions:      (lat: 69, lon: 248)
Coordinates:
  * lat          (lat) float64 552B 51.88 52.12 52.38 ... 68.38 68.62 68.88
  * lon          (lon) float64 2kB -66.88 -66.62 -66.38 ... -5.625 -5.375 -5.125
Data variables:
    probability  (lat, lon) float64 137kB nan nan nan nan ... nan nan nan nan
[14]:
# Plotting the Lagrangian probability in discrete 2-D longitude-latitude bins:
traj_geo.summary_data.probability.plot()

[14]:
<matplotlib.collections.QuadMesh at 0x7f283d0e3df0>
../_images/tutorials_tutorial_2_advanced_analysis_with_TrajFrames_22_1.png