Langmuir turbulence example

This example implements a Langmuir turbulence simulation similar to the one reported in section 4 of

Wagner et al., "Near-inertial waves and turbulence driven by the growth of swell", Journal of Physical Oceanography (2021)

This example demonstrates

  • How to run large eddy simulations with surface wave effects via the Craik-Leibovich approximation.

  • How to specify time- and horizontally-averaged output.

Install dependencies

First let's make sure we have all required packages installed.

using Pkg
pkg"add Oceananigans, CairoMakie, CUDA"
using Oceananigans
using Oceananigans.Units: minute, minutes, hours
using CUDA

Model set-up

To build the model, we specify the grid, Stokes drift, boundary conditions, and Coriolis parameter.

Domain and numerical grid specification

We use a modest resolution and the same total extent as Wagner et al. (2021),

grid = RectilinearGrid(GPU(), size=(128, 128, 64), extent=(128, 128, 64))
128×128×64 RectilinearGrid{Float64, Periodic, Periodic, Bounded} on CUDAGPU with 3×3×3 halo
├── Periodic x ∈ [0.0, 128.0) regularly spaced with Δx=1.0
├── Periodic y ∈ [0.0, 128.0) regularly spaced with Δy=1.0
└── Bounded  z ∈ [-64.0, 0.0] regularly spaced with Δz=1.0

The Stokes Drift profile

The surface wave Stokes drift profile prescribed in Wagner et al. (2021), corresponds to a 'monochromatic' (that is, single-frequency) wave field.

A monochromatic wave field is characterized by its wavelength and amplitude (half the distance from wave crest to wave trough), which determine the wave frequency and the vertical scale of the Stokes drift profile.

g = Oceananigans.defaults.gravitational_acceleration

amplitude = 0.8 # m
wavelength = 60  # m
wavenumber = 2π / wavelength # m⁻¹
frequency = sqrt(g * wavenumber) # s⁻¹

# The vertical scale over which the Stokes drift of a monochromatic surface wave
# decays away from the surface is `1/2wavenumber`, or
const vertical_scale = wavelength / 4π

# Stokes drift velocity at the surface
const Uˢ = amplitude^2 * wavenumber * frequency # m s⁻¹
0.06791774197745354

The const declarations ensure that Stokes drift functions compile on the GPU. To run this example on the CPU, replace GPU() with CPU() in the RectilinearGrid constructor above.

The Stokes drift profile is

uˢ(z) = Uˢ * exp(z / vertical_scale)
uˢ (generic function with 1 method)

and its z-derivative is

∂z_uˢ(z, t) = 1 / vertical_scale * Uˢ * exp(z / vertical_scale)
∂z_uˢ (generic function with 1 method)
The Craik-Leibovich equations in Oceananigans

Oceananigans implements the Craik-Leibovich approximation for surface wave effects using the Lagrangian-mean velocity field as its prognostic momentum variable. In other words, model.velocities.u is the Lagrangian-mean $x$-velocity beneath surface waves. This differs from models that use the Eulerian-mean velocity field as a prognostic variable, but has the advantage that $u$ accounts for the total advection of tracers and momentum, and that $u = v = w = 0$ is a steady solution even when Coriolis forces are present. See the physics documentation for more information.

Finally, we note that the time-derivative of the Stokes drift must be provided if the Stokes drift and surface wave field undergoes forced changes in time. In this example, the Stokes drift is constant and thus the time-derivative of the Stokes drift is 0.

Boundary conditions

At the surface $z = 0$, Wagner et al. (2021) impose

τx = -3.72e-5 # m² s⁻², surface kinematic momentum flux
u_boundary_conditions = FieldBoundaryConditions(top = FluxBoundaryCondition(τx))
Oceananigans.FieldBoundaryConditions, with boundary conditions
├── west: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── east: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── south: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── north: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── bottom: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── top: FluxBoundaryCondition: -3.72e-5
└── immersed: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)

Wagner et al. (2021) impose a linear buoyancy gradient at the bottom along with a weak, destabilizing flux of buoyancy at the surface to faciliate spin-up from rest.

Jᵇ = 2.307e-8 # m² s⁻³, surface buoyancy flux
N² = 1.936e-5 # s⁻², initial and bottom buoyancy gradient

b_boundary_conditions = FieldBoundaryConditions(top = FluxBoundaryCondition(Jᵇ),
                                                bottom = GradientBoundaryCondition(N²))
Oceananigans.FieldBoundaryConditions, with boundary conditions
├── west: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── east: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── south: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── north: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── bottom: GradientBoundaryCondition: 1.936e-5
├── top: FluxBoundaryCondition: 2.307e-8
└── immersed: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
The flux convention in Oceananigans

Note that Oceananigans uses "positive upward" conventions for all fluxes. In consequence, a negative flux at the surface drives positive velocities, and a positive flux of buoyancy drives cooling.

Coriolis parameter

Wagner et al. (2021) use

coriolis = FPlane(f=1e-4) # s⁻¹
FPlane{Float64}(f=0.0001)

which is typical for mid-latitudes on Earth.

Model instantiation

We are ready to build the model. We use a fifth-order Weighted Essentially Non-Oscillatory (WENO) advection scheme and the AnisotropicMinimumDissipation model for large eddy simulation. Because our Stokes drift does not vary in $x, y$, we use UniformStokesDrift, which expects Stokes drift functions of $z, t$ only.

model = NonhydrostaticModel(grid; coriolis,
                            advection = WENO(order=9),
                            tracers = :b,
                            buoyancy = BuoyancyTracer(),
                            stokes_drift = UniformStokesDrift(∂z_uˢ=∂z_uˢ),
                            boundary_conditions = (u=u_boundary_conditions, b=b_boundary_conditions))
NonhydrostaticModel{CUDAGPU, RectilinearGrid}(time = 0 seconds, iteration = 0)
├── grid: 128×128×64 RectilinearGrid{Float64, Periodic, Periodic, Bounded} on CUDAGPU with 5×5×5 halo
├── timestepper: RungeKutta3TimeStepper
├── advection scheme: WENO{5, Float64, Float32}(order=9)
├── tracers: b
├── closure: Nothing
├── buoyancy: BuoyancyTracer with ĝ = NegativeZDirection()
└── coriolis: FPlane{Float64}(f=0.0001)

Initial conditions

We make use of random noise concentrated in the upper 4 meters for buoyancy and velocity initial conditions,

Ξ(z) = randn() * exp(z / 4)

Our initial condition for buoyancy consists of a surface mixed layer 33 m deep, a deep linear stratification, plus noise,

initial_mixed_layer_depth = 33 # m
stratification(z) = z < - initial_mixed_layer_depth ? N² * z : N² * (-initial_mixed_layer_depth)

bᵢ(x, y, z) = stratification(z) + 1e-1 * Ξ(z) * N² * model.grid.Lz
bᵢ (generic function with 1 method)

The simulation we reproduce from Wagner et al. (2021) is zero Lagrangian-mean velocity. This initial condition is consistent with a wavy, quiescent ocean suddenly impacted by winds. To this quiescent state we add noise scaled by the friction velocity to $u$ and $w$.

u★ = sqrt(abs(τx))
uᵢ(x, y, z) = u★ * 1e-1 * Ξ(z)
wᵢ(x, y, z) = u★ * 1e-1 * Ξ(z)

set!(model, u=uᵢ, w=wᵢ, b=bᵢ)

Setting up the simulation

simulation = Simulation(model, Δt=45.0, stop_time=4hours)
Simulation of NonhydrostaticModel{CUDAGPU, RectilinearGrid}(time = 0 seconds, iteration = 0)
├── Next time step: 45 seconds
├── run_wall_time: 0 seconds
├── run_wall_time / iteration: NaN days
├── stop_time: 4 hours
├── stop_iteration: Inf
├── wall_time_limit: Inf
├── minimum_relative_step: 0.0
├── callbacks: OrderedDict with 4 entries:
│   ├── stop_time_exceeded => Callback of stop_time_exceeded on IterationInterval(1)
│   ├── stop_iteration_exceeded => Callback of stop_iteration_exceeded on IterationInterval(1)
│   ├── wall_time_limit_exceeded => Callback of wall_time_limit_exceeded on IterationInterval(1)
│   └── nan_checker => Callback of NaNChecker for u on IterationInterval(100)
└── output_writers: OrderedDict with no entries

We use the TimeStepWizard for adaptive time-stepping with a Courant-Freidrichs-Lewy (CFL) number of 1.0,

conjure_time_step_wizard!(simulation, cfl=1.0, max_Δt=1minute)

Nice progress messaging

We define a function that prints a helpful message with maximum absolute value of $u, v, w$ and the current wall clock time.

using Printf

function progress(simulation)
    u, v, w = simulation.model.velocities

    # Print a progress message
    msg = @sprintf("i: %04d, t: %s, Δt: %s, umax = (%.1e, %.1e, %.1e) ms⁻¹, wall time: %s\n",
                   iteration(simulation),
                   prettytime(time(simulation)),
                   prettytime(simulation.Δt),
                   maximum(abs, u), maximum(abs, v), maximum(abs, w),
                   prettytime(simulation.run_wall_time))

    @info msg

    return nothing
end

simulation.callbacks[:progress] = Callback(progress, IterationInterval(20))
Callback of progress on IterationInterval(20)

Output

A field writer

We set up an output writer for the simulation that saves all velocity fields, tracer fields, and the subgrid turbulent diffusivity.

output_interval = 5minutes

fields_to_output = merge(model.velocities, model.tracers)

simulation.output_writers[:fields] =
    JLD2Writer(model, fields_to_output,
               schedule = TimeInterval(output_interval),
               filename = "langmuir_turbulence_fields.jld2",
               overwrite_existing = true)
JLD2Writer scheduled on TimeInterval(5 minutes):
├── filepath: langmuir_turbulence_fields.jld2
├── 4 outputs: (u, v, w, b)
├── array_type: Array{Float32}
├── including: [:grid, :coriolis, :buoyancy, :closure]
├── file_splitting: NoFileSplitting
└── file size: 0 bytes (file not yet created)

An "averages" writer

We also set up output of time- and horizontally-averaged velocity field and momentum fluxes.

u, v, w = model.velocities
b = model.tracers.b

 U = Average(u, dims=(1, 2))
 V = Average(v, dims=(1, 2))
 B = Average(b, dims=(1, 2))
wu = Average(w * u, dims=(1, 2))
wv = Average(w * v, dims=(1, 2))

simulation.output_writers[:averages] =
    JLD2Writer(model, (; U, V, B, wu, wv),
               schedule = AveragedTimeInterval(output_interval, window=2minutes),
               filename = "langmuir_turbulence_averages.jld2",
               overwrite_existing = true)
JLD2Writer scheduled on TimeInterval(5 minutes):
├── filepath: langmuir_turbulence_averages.jld2
├── 5 outputs: (U, V, B, wu, wv) averaged on AveragedTimeInterval(window=2 minutes, stride=1, interval=5 minutes)
├── array_type: Array{Float32}
├── including: [:grid, :coriolis, :buoyancy, :closure]
├── file_splitting: NoFileSplitting
└── file size: 0 bytes (file not yet created)

Running the simulation

This part is easy,

run!(simulation)
[ Info: Initializing simulation...
[ Info: i: 0000, t: 0 seconds, Δt: 49.500 seconds, umax = (1.6e-03, 8.0e-04, 1.4e-03) ms⁻¹, wall time: 0 seconds
[ Info:     ... simulation initialization complete (12.115 seconds)
[ Info: Executing initial time step...
[ Info:     ... initial time step complete (2.888 seconds).
[ Info: i: 0020, t: 11.238 minutes, Δt: 19.473 seconds, umax = (3.5e-02, 1.3e-02, 2.2e-02) ms⁻¹, wall time: 16.201 seconds
[ Info: i: 0040, t: 16.750 minutes, Δt: 12.623 seconds, umax = (5.1e-02, 2.1e-02, 2.4e-02) ms⁻¹, wall time: 16.626 seconds
[ Info: i: 0060, t: 20.777 minutes, Δt: 10.872 seconds, umax = (6.3e-02, 2.8e-02, 3.2e-02) ms⁻¹, wall time: 17.143 seconds
[ Info: i: 0080, t: 24.465 minutes, Δt: 10.762 seconds, umax = (6.4e-02, 3.2e-02, 3.3e-02) ms⁻¹, wall time: 17.504 seconds
[ Info: i: 0100, t: 28.031 minutes, Δt: 10.583 seconds, umax = (6.2e-02, 3.2e-02, 3.4e-02) ms⁻¹, wall time: 17.940 seconds
[ Info: i: 0120, t: 31.457 minutes, Δt: 10.694 seconds, umax = (6.1e-02, 3.3e-02, 2.8e-02) ms⁻¹, wall time: 18.382 seconds
[ Info: i: 0140, t: 35 minutes, Δt: 10.253 seconds, umax = (6.6e-02, 3.3e-02, 3.0e-02) ms⁻¹, wall time: 18.798 seconds
[ Info: i: 0160, t: 38.474 minutes, Δt: 10.487 seconds, umax = (7.1e-02, 3.7e-02, 3.7e-02) ms⁻¹, wall time: 19.236 seconds
[ Info: i: 0180, t: 41.773 minutes, Δt: 9.997 seconds, umax = (7.2e-02, 3.7e-02, 3.6e-02) ms⁻¹, wall time: 19.674 seconds
[ Info: i: 0200, t: 44.882 minutes, Δt: 9.294 seconds, umax = (7.4e-02, 3.7e-02, 3.5e-02) ms⁻¹, wall time: 20.089 seconds
[ Info: i: 0220, t: 47.815 minutes, Δt: 9.101 seconds, umax = (7.5e-02, 3.8e-02, 3.6e-02) ms⁻¹, wall time: 20.518 seconds
[ Info: i: 0240, t: 50.723 minutes, Δt: 9.335 seconds, umax = (7.3e-02, 4.2e-02, 3.7e-02) ms⁻¹, wall time: 21.000 seconds
[ Info: i: 0260, t: 53.787 minutes, Δt: 8.939 seconds, umax = (7.9e-02, 4.4e-02, 3.6e-02) ms⁻¹, wall time: 21.381 seconds
[ Info: i: 0280, t: 56.671 minutes, Δt: 8.707 seconds, umax = (8.3e-02, 4.1e-02, 3.6e-02) ms⁻¹, wall time: 21.815 seconds
[ Info: i: 0300, t: 59.593 minutes, Δt: 8.146 seconds, umax = (8.0e-02, 4.4e-02, 3.8e-02) ms⁻¹, wall time: 22.234 seconds
[ Info: i: 0320, t: 1.039 hours, Δt: 8.650 seconds, umax = (8.0e-02, 5.0e-02, 3.4e-02) ms⁻¹, wall time: 22.669 seconds
[ Info: i: 0340, t: 1.086 hours, Δt: 7.389 seconds, umax = (8.4e-02, 4.8e-02, 3.8e-02) ms⁻¹, wall time: 23.236 seconds
[ Info: i: 0360, t: 1.128 hours, Δt: 7.872 seconds, umax = (8.1e-02, 4.4e-02, 3.6e-02) ms⁻¹, wall time: 23.535 seconds
[ Info: i: 0380, t: 1.171 hours, Δt: 8.196 seconds, umax = (8.3e-02, 4.7e-02, 3.8e-02) ms⁻¹, wall time: 24.086 seconds
[ Info: i: 0400, t: 1.217 hours, Δt: 8.493 seconds, umax = (8.9e-02, 4.7e-02, 4.0e-02) ms⁻¹, wall time: 24.406 seconds
[ Info: i: 0420, t: 1.261 hours, Δt: 7.583 seconds, umax = (9.0e-02, 4.7e-02, 3.8e-02) ms⁻¹, wall time: 24.911 seconds
[ Info: i: 0440, t: 1.303 hours, Δt: 8.061 seconds, umax = (8.7e-02, 4.6e-02, 4.1e-02) ms⁻¹, wall time: 25.287 seconds
[ Info: i: 0460, t: 1.346 hours, Δt: 7.530 seconds, umax = (9.2e-02, 5.3e-02, 3.8e-02) ms⁻¹, wall time: 25.766 seconds
[ Info: i: 0480, t: 1.388 hours, Δt: 7.352 seconds, umax = (9.2e-02, 5.3e-02, 3.8e-02) ms⁻¹, wall time: 26.170 seconds
[ Info: i: 0500, t: 1.430 hours, Δt: 7.465 seconds, umax = (9.1e-02, 5.1e-02, 4.2e-02) ms⁻¹, wall time: 26.652 seconds
[ Info: i: 0520, t: 1.470 hours, Δt: 7.122 seconds, umax = (9.4e-02, 6.0e-02, 4.4e-02) ms⁻¹, wall time: 27.058 seconds
[ Info: i: 0540, t: 1.508 hours, Δt: 6.979 seconds, umax = (9.1e-02, 5.9e-02, 4.6e-02) ms⁻¹, wall time: 27.585 seconds
[ Info: i: 0560, t: 1.548 hours, Δt: 7.400 seconds, umax = (9.3e-02, 5.4e-02, 4.0e-02) ms⁻¹, wall time: 27.952 seconds
[ Info: i: 0580, t: 1.587 hours, Δt: 7.298 seconds, umax = (9.8e-02, 6.0e-02, 4.6e-02) ms⁻¹, wall time: 28.514 seconds
[ Info: i: 0600, t: 1.628 hours, Δt: 6.997 seconds, umax = (9.7e-02, 5.1e-02, 4.2e-02) ms⁻¹, wall time: 28.845 seconds
[ Info: i: 0620, t: 1.667 hours, Δt: 6.682 seconds, umax = (9.6e-02, 5.5e-02, 4.3e-02) ms⁻¹, wall time: 29.287 seconds
[ Info: i: 0640, t: 1.706 hours, Δt: 7.165 seconds, umax = (1.0e-01, 5.5e-02, 5.0e-02) ms⁻¹, wall time: 29.753 seconds
[ Info: i: 0660, t: 1.745 hours, Δt: 7.016 seconds, umax = (9.2e-02, 5.5e-02, 4.4e-02) ms⁻¹, wall time: 30.200 seconds
[ Info: i: 0680, t: 1.783 hours, Δt: 7.356 seconds, umax = (9.5e-02, 5.6e-02, 4.3e-02) ms⁻¹, wall time: 30.671 seconds
[ Info: i: 0700, t: 1.822 hours, Δt: 6.632 seconds, umax = (9.7e-02, 6.4e-02, 4.2e-02) ms⁻¹, wall time: 31.121 seconds
[ Info: i: 0720, t: 1.860 hours, Δt: 6.832 seconds, umax = (1.1e-01, 5.9e-02, 4.2e-02) ms⁻¹, wall time: 31.598 seconds
[ Info: i: 0740, t: 1.897 hours, Δt: 6.466 seconds, umax = (1.0e-01, 6.5e-02, 4.2e-02) ms⁻¹, wall time: 32.046 seconds
[ Info: i: 0760, t: 1.933 hours, Δt: 6.519 seconds, umax = (9.8e-02, 6.1e-02, 4.5e-02) ms⁻¹, wall time: 32.535 seconds
[ Info: i: 0780, t: 1.970 hours, Δt: 7.006 seconds, umax = (1.1e-01, 6.2e-02, 4.7e-02) ms⁻¹, wall time: 32.970 seconds
[ Info: i: 0800, t: 2.007 hours, Δt: 6.874 seconds, umax = (1.0e-01, 6.1e-02, 4.4e-02) ms⁻¹, wall time: 33.520 seconds
[ Info: i: 0820, t: 2.045 hours, Δt: 6.743 seconds, umax = (1.0e-01, 6.0e-02, 4.3e-02) ms⁻¹, wall time: 33.902 seconds
[ Info: i: 0840, t: 2.082 hours, Δt: 6.486 seconds, umax = (1.0e-01, 5.9e-02, 4.4e-02) ms⁻¹, wall time: 34.360 seconds
[ Info: i: 0860, t: 2.119 hours, Δt: 6.354 seconds, umax = (1.0e-01, 6.1e-02, 4.5e-02) ms⁻¹, wall time: 34.829 seconds
[ Info: i: 0880, t: 2.154 hours, Δt: 6.296 seconds, umax = (1.0e-01, 6.5e-02, 4.2e-02) ms⁻¹, wall time: 35.286 seconds
[ Info: i: 0900, t: 2.188 hours, Δt: 6.566 seconds, umax = (1.0e-01, 6.7e-02, 4.7e-02) ms⁻¹, wall time: 35.768 seconds
[ Info: i: 0920, t: 2.224 hours, Δt: 6.174 seconds, umax = (1.0e-01, 5.8e-02, 4.5e-02) ms⁻¹, wall time: 36.221 seconds
[ Info: i: 0940, t: 2.259 hours, Δt: 6.369 seconds, umax = (1.0e-01, 5.9e-02, 4.3e-02) ms⁻¹, wall time: 36.784 seconds
[ Info: i: 0960, t: 2.294 hours, Δt: 6.256 seconds, umax = (1.0e-01, 6.2e-02, 4.3e-02) ms⁻¹, wall time: 37.187 seconds
[ Info: i: 0980, t: 2.329 hours, Δt: 6.366 seconds, umax = (1.0e-01, 6.4e-02, 4.3e-02) ms⁻¹, wall time: 37.648 seconds
[ Info: i: 1000, t: 2.365 hours, Δt: 6.262 seconds, umax = (1.0e-01, 6.2e-02, 3.9e-02) ms⁻¹, wall time: 38.127 seconds
[ Info: i: 1020, t: 2.400 hours, Δt: 6.648 seconds, umax = (1.0e-01, 6.2e-02, 4.7e-02) ms⁻¹, wall time: 38.585 seconds
[ Info: i: 1040, t: 2.437 hours, Δt: 6.681 seconds, umax = (1.0e-01, 6.7e-02, 4.9e-02) ms⁻¹, wall time: 39.065 seconds
[ Info: i: 1060, t: 2.473 hours, Δt: 6.276 seconds, umax = (1.0e-01, 7.2e-02, 4.5e-02) ms⁻¹, wall time: 39.516 seconds
[ Info: i: 1080, t: 2.507 hours, Δt: 5.833 seconds, umax = (1.1e-01, 7.7e-02, 5.0e-02) ms⁻¹, wall time: 40.073 seconds
[ Info: i: 1100, t: 2.540 hours, Δt: 6.419 seconds, umax = (1.1e-01, 6.6e-02, 5.0e-02) ms⁻¹, wall time: 40.459 seconds
[ Info: i: 1120, t: 2.575 hours, Δt: 6.063 seconds, umax = (1.1e-01, 7.3e-02, 4.9e-02) ms⁻¹, wall time: 40.921 seconds
[ Info: i: 1140, t: 2.608 hours, Δt: 6.148 seconds, umax = (1.0e-01, 6.8e-02, 4.5e-02) ms⁻¹, wall time: 41.416 seconds
[ Info: i: 1160, t: 2.642 hours, Δt: 6.131 seconds, umax = (1.1e-01, 7.0e-02, 4.5e-02) ms⁻¹, wall time: 41.879 seconds
[ Info: i: 1180, t: 2.675 hours, Δt: 6.240 seconds, umax = (1.1e-01, 6.7e-02, 4.7e-02) ms⁻¹, wall time: 42.405 seconds
[ Info: i: 1200, t: 2.711 hours, Δt: 6.098 seconds, umax = (1.1e-01, 6.5e-02, 4.5e-02) ms⁻¹, wall time: 42.806 seconds
[ Info: i: 1220, t: 2.744 hours, Δt: 6.337 seconds, umax = (1.1e-01, 7.0e-02, 4.9e-02) ms⁻¹, wall time: 43.271 seconds
[ Info: i: 1240, t: 2.778 hours, Δt: 6.254 seconds, umax = (1.1e-01, 6.7e-02, 4.7e-02) ms⁻¹, wall time: 43.764 seconds
[ Info: i: 1260, t: 2.813 hours, Δt: 6.635 seconds, umax = (1.1e-01, 6.4e-02, 4.2e-02) ms⁻¹, wall time: 44.225 seconds
[ Info: i: 1280, t: 2.847 hours, Δt: 6.021 seconds, umax = (1.1e-01, 6.9e-02, 4.3e-02) ms⁻¹, wall time: 44.709 seconds
[ Info: i: 1300, t: 2.880 hours, Δt: 6.252 seconds, umax = (1.1e-01, 7.1e-02, 4.2e-02) ms⁻¹, wall time: 45.153 seconds
[ Info: i: 1320, t: 2.915 hours, Δt: 6.031 seconds, umax = (1.0e-01, 7.8e-02, 4.5e-02) ms⁻¹, wall time: 45.616 seconds
[ Info: i: 1340, t: 2.947 hours, Δt: 6.212 seconds, umax = (1.0e-01, 7.4e-02, 4.8e-02) ms⁻¹, wall time: 46.108 seconds
[ Info: i: 1360, t: 2.982 hours, Δt: 6.684 seconds, umax = (1.0e-01, 6.9e-02, 4.8e-02) ms⁻¹, wall time: 46.568 seconds
[ Info: i: 1380, t: 3.018 hours, Δt: 6.012 seconds, umax = (1.1e-01, 7.6e-02, 5.3e-02) ms⁻¹, wall time: 47.057 seconds
[ Info: i: 1400, t: 3.051 hours, Δt: 6.043 seconds, umax = (1.1e-01, 7.1e-02, 4.6e-02) ms⁻¹, wall time: 47.504 seconds
[ Info: i: 1420, t: 3.083 hours, Δt: 6.193 seconds, umax = (1.1e-01, 7.5e-02, 4.6e-02) ms⁻¹, wall time: 47.971 seconds
[ Info: i: 1440, t: 3.118 hours, Δt: 5.901 seconds, umax = (1.0e-01, 7.6e-02, 4.8e-02) ms⁻¹, wall time: 48.456 seconds
[ Info: i: 1460, t: 3.151 hours, Δt: 6.418 seconds, umax = (1.1e-01, 7.4e-02, 4.6e-02) ms⁻¹, wall time: 48.920 seconds
[ Info: i: 1480, t: 3.186 hours, Δt: 6.173 seconds, umax = (1.2e-01, 6.9e-02, 4.9e-02) ms⁻¹, wall time: 49.402 seconds
[ Info: i: 1500, t: 3.220 hours, Δt: 6.078 seconds, umax = (1.2e-01, 7.3e-02, 5.2e-02) ms⁻¹, wall time: 49.848 seconds
[ Info: i: 1520, t: 3.253 hours, Δt: 6.328 seconds, umax = (1.1e-01, 6.8e-02, 4.6e-02) ms⁻¹, wall time: 50.452 seconds
[ Info: i: 1540, t: 3.288 hours, Δt: 6.194 seconds, umax = (1.1e-01, 7.0e-02, 4.8e-02) ms⁻¹, wall time: 50.798 seconds
[ Info: i: 1560, t: 3.322 hours, Δt: 6.064 seconds, umax = (1.1e-01, 7.4e-02, 4.8e-02) ms⁻¹, wall time: 51.261 seconds
[ Info: i: 1580, t: 3.354 hours, Δt: 5.958 seconds, umax = (1.2e-01, 7.5e-02, 4.7e-02) ms⁻¹, wall time: 51.761 seconds
[ Info: i: 1600, t: 3.387 hours, Δt: 5.170 seconds, umax = (1.1e-01, 7.5e-02, 4.4e-02) ms⁻¹, wall time: 52.213 seconds
[ Info: i: 1620, t: 3.416 hours, Δt: 5.720 seconds, umax = (1.1e-01, 7.0e-02, 4.7e-02) ms⁻¹, wall time: 52.673 seconds
[ Info: i: 1640, t: 3.447 hours, Δt: 5.570 seconds, umax = (1.1e-01, 7.4e-02, 5.1e-02) ms⁻¹, wall time: 53.163 seconds
[ Info: i: 1660, t: 3.478 hours, Δt: 5.642 seconds, umax = (1.1e-01, 7.3e-02, 4.7e-02) ms⁻¹, wall time: 53.625 seconds
[ Info: i: 1680, t: 3.509 hours, Δt: 5.073 seconds, umax = (1.1e-01, 8.1e-02, 4.5e-02) ms⁻¹, wall time: 54.137 seconds
[ Info: i: 1700, t: 3.539 hours, Δt: 6.003 seconds, umax = (1.1e-01, 7.2e-02, 5.3e-02) ms⁻¹, wall time: 54.563 seconds
[ Info: i: 1720, t: 3.571 hours, Δt: 6.210 seconds, umax = (1.0e-01, 7.6e-02, 5.1e-02) ms⁻¹, wall time: 55.027 seconds
[ Info: i: 1740, t: 3.604 hours, Δt: 5.778 seconds, umax = (1.1e-01, 7.8e-02, 4.9e-02) ms⁻¹, wall time: 55.508 seconds
[ Info: i: 1760, t: 3.637 hours, Δt: 5.980 seconds, umax = (1.1e-01, 7.1e-02, 5.0e-02) ms⁻¹, wall time: 55.955 seconds
[ Info: i: 1780, t: 3.668 hours, Δt: 5.587 seconds, umax = (1.1e-01, 8.0e-02, 5.4e-02) ms⁻¹, wall time: 56.575 seconds
[ Info: i: 1800, t: 3.701 hours, Δt: 5.038 seconds, umax = (1.2e-01, 8.0e-02, 4.7e-02) ms⁻¹, wall time: 56.902 seconds
[ Info: i: 1820, t: 3.728 hours, Δt: 5.523 seconds, umax = (1.3e-01, 7.6e-02, 4.5e-02) ms⁻¹, wall time: 57.363 seconds
[ Info: i: 1840, t: 3.760 hours, Δt: 6.071 seconds, umax = (1.1e-01, 8.1e-02, 4.5e-02) ms⁻¹, wall time: 57.863 seconds
[ Info: i: 1860, t: 3.792 hours, Δt: 5.567 seconds, umax = (1.1e-01, 8.1e-02, 4.7e-02) ms⁻¹, wall time: 58.291 seconds
[ Info: i: 1880, t: 3.823 hours, Δt: 5.853 seconds, umax = (1.2e-01, 8.2e-02, 4.9e-02) ms⁻¹, wall time: 58.755 seconds
[ Info: i: 1900, t: 3.854 hours, Δt: 5.858 seconds, umax = (1.0e-01, 9.3e-02, 5.0e-02) ms⁻¹, wall time: 59.237 seconds
[ Info: i: 1920, t: 3.886 hours, Δt: 5.855 seconds, umax = (1.2e-01, 8.0e-02, 5.6e-02) ms⁻¹, wall time: 59.683 seconds
[ Info: i: 1940, t: 3.918 hours, Δt: 6.192 seconds, umax = (1.1e-01, 7.7e-02, 5.9e-02) ms⁻¹, wall time: 1.005 minutes
[ Info: i: 1960, t: 3.952 hours, Δt: 5.486 seconds, umax = (1.1e-01, 8.3e-02, 6.0e-02) ms⁻¹, wall time: 1.010 minutes
[ Info: i: 1980, t: 3.981 hours, Δt: 5.429 seconds, umax = (1.1e-01, 7.8e-02, 5.6e-02) ms⁻¹, wall time: 1.018 minutes
[ Info: Simulation is stopping after running for 1.024 minutes.
[ Info: Simulation time 4 hours equals or exceeds stop time 4 hours.

Making a neat movie

We look at the results by loading data from file with FieldTimeSeries, and plotting vertical slices of $u$ and $w$, and a horizontal slice of $w$ to look for Langmuir cells.

using CairoMakie

time_series = (;
     w = FieldTimeSeries("langmuir_turbulence_fields.jld2", "w"),
     u = FieldTimeSeries("langmuir_turbulence_fields.jld2", "u"),
     B = FieldTimeSeries("langmuir_turbulence_averages.jld2", "B"),
     U = FieldTimeSeries("langmuir_turbulence_averages.jld2", "U"),
     V = FieldTimeSeries("langmuir_turbulence_averages.jld2", "V"),
    wu = FieldTimeSeries("langmuir_turbulence_averages.jld2", "wu"),
    wv = FieldTimeSeries("langmuir_turbulence_averages.jld2", "wv"))

times = time_series.w.times

We are now ready to animate using Makie. We use Makie's Observable to animate the data. To dive into how Observables work we refer to Makie.jl's Documentation.

n = Observable(1)

wxy_title = @lift string("w(x, y, t) at z=-8 m and t = ", prettytime(times[$n]))
wxz_title = @lift string("w(x, z, t) at y=0 m and t = ", prettytime(times[$n]))
uxz_title = @lift string("u(x, z, t) at y=0 m and t = ", prettytime(times[$n]))

fig = Figure(size = (850, 850))

ax_B = Axis(fig[1, 4];
            xlabel = "Buoyancy (m s⁻²)",
            ylabel = "z (m)")

ax_U = Axis(fig[2, 4];
            xlabel = "Velocities (m s⁻¹)",
            ylabel = "z (m)",
            limits = ((-0.07, 0.07), nothing))

ax_fluxes = Axis(fig[3, 4];
                 xlabel = "Momentum fluxes (m² s⁻²)",
                 ylabel = "z (m)",
                 limits = ((-3.5e-5, 3.5e-5), nothing))

ax_wxy = Axis(fig[1, 1:2];
              xlabel = "x (m)",
              ylabel = "y (m)",
              aspect = DataAspect(),
              limits = ((0, grid.Lx), (0, grid.Ly)),
              title = wxy_title)

ax_wxz = Axis(fig[2, 1:2];
              xlabel = "x (m)",
              ylabel = "z (m)",
              aspect = AxisAspect(2),
              limits = ((0, grid.Lx), (-grid.Lz, 0)),
              title = wxz_title)

ax_uxz = Axis(fig[3, 1:2];
              xlabel = "x (m)",
              ylabel = "z (m)",
              aspect = AxisAspect(2),
              limits = ((0, grid.Lx), (-grid.Lz, 0)),
              title = uxz_title)


wₙ = @lift time_series.w[$n]
uₙ = @lift time_series.u[$n]
Bₙ = @lift view(time_series.B[$n], 1, 1, :)
Uₙ = @lift view(time_series.U[$n], 1, 1, :)
Vₙ = @lift view(time_series.V[$n], 1, 1, :)
wuₙ = @lift view(time_series.wu[$n], 1, 1, :)
wvₙ = @lift view(time_series.wv[$n], 1, 1, :)

k = searchsortedfirst(znodes(grid, Face(); with_halos=true), -8)
wxyₙ = @lift view(time_series.w[$n], :, :, k)
wxzₙ = @lift view(time_series.w[$n], :, 1, :)
uxzₙ = @lift view(time_series.u[$n], :, 1, :)

wlims = (-0.03, 0.03)
ulims = (-0.05, 0.05)

lines!(ax_B, Bₙ)

lines!(ax_U, Uₙ; label = L"\bar{u}")
lines!(ax_U, Vₙ; label = L"\bar{v}")
axislegend(ax_U; position = :rb)

lines!(ax_fluxes, wuₙ; label = L"mean $wu$")
lines!(ax_fluxes, wvₙ; label = L"mean $wv$")
axislegend(ax_fluxes; position = :rb)

hm_wxy = heatmap!(ax_wxy, wxyₙ;
                  colorrange = wlims,
                  colormap = :balance)

Colorbar(fig[1, 3], hm_wxy; label = "m s⁻¹")

hm_wxz = heatmap!(ax_wxz, wxzₙ;
                  colorrange = wlims,
                  colormap = :balance)

Colorbar(fig[2, 3], hm_wxz; label = "m s⁻¹")

ax_uxz = heatmap!(ax_uxz, uxzₙ;
                  colorrange = ulims,
                  colormap = :balance)

Colorbar(fig[3, 3], ax_uxz; label = "m s⁻¹")

fig

And, finally, we record a movie.

frames = 1:length(times)

CairoMakie.record(fig, "langmuir_turbulence.mp4", frames, framerate=8) do i
    n[] = i
end


Julia version and environment information

This example was executed with the following version of Julia:

using InteractiveUtils: versioninfo
versioninfo()
Julia Version 1.12.4
Commit 01a2eadb047 (2026-01-06 16:56 UTC)
Build Info:
  Official https://julialang.org release
Platform Info:
  OS: Linux (x86_64-linux-gnu)
  CPU: 128 × AMD EPYC 9374F 32-Core Processor
  WORD_SIZE: 64
  LLVM: libLLVM-18.1.7 (ORCJIT, znver4)
  GC: Built with stock GC
Threads: 1 default, 1 interactive, 1 GC (on 128 virtual cores)
Environment:
  LD_LIBRARY_PATH = 
  JULIA_PKG_SERVER_REGISTRY_PREFERENCE = eager
  JULIA_DEPOT_PATH = /var/lib/buildkite-agent/.julia-oceananigans
  JULIA_PROJECT = /var/lib/buildkite-agent/Oceananigans.jl-28584/docs/
  JULIA_VERSION = 1.12.4
  JULIA_LOAD_PATH = @:@v#.#:@stdlib
  JULIA_VERSION_ENZYME = 1.10.10
  JULIA_PYTHONCALL_EXE = /var/lib/buildkite-agent/Oceananigans.jl-28584/docs/.CondaPkg/.pixi/envs/default/bin/python
  JULIA_DEBUG = Literate

These were the top-level packages installed in the environment:

import Pkg
Pkg.status()
Status `~/Oceananigans.jl-28584/docs/Project.toml`
  [79e6a3ab] Adapt v4.4.0
  [052768ef] CUDA v5.9.6
  [13f3f980] CairoMakie v0.15.8
  [e30172f5] Documenter v1.16.1
  [daee34ce] DocumenterCitations v1.4.1
  [033835bb] JLD2 v0.6.3
  [63c18a36] KernelAbstractions v0.9.39
  [98b081ad] Literate v2.21.0
  [da04e1cc] MPI v0.20.23
  [85f8d34a] NCDatasets v0.14.10
  [9e8cae18] Oceananigans v0.104.0 `..`
  [f27b6e38] Polynomials v4.1.0
  [6038ab10] Rotations v1.7.1
  [d496a93d] SeawaterPolynomials v0.3.10
  [09ab397b] StructArrays v0.7.2
  [bdfc003b] TimesDates v0.3.3
  [2e0b0046] XESMF v0.1.6
  [b77e0a4c] InteractiveUtils v1.11.0
  [37e2e46d] LinearAlgebra v1.12.0
  [44cfe95a] Pkg v1.12.1

This page was generated using Literate.jl.