Skip to content

Wind- and convection-driven mixing in an ocean surface boundary layer

This example simulates mixing by three-dimensional turbulence in an ocean surface boundary layer driven by atmospheric winds and convection. It demonstrates:

  • How to set-up a grid with varying spacing in the vertical direction

  • How to use the SeawaterBuoyancy model for buoyancy with TEOS10EquationOfState.

  • How to use a turbulence closure for large eddy simulation.

  • How to use a function to impose a boundary condition.

Install dependencies

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

julia
using Pkg
pkg"add Oceananigans, CairoMakie, SeawaterPolynomials, CUDA"

We start by importing all of the packages and functions that we'll need for this example.

julia
using Oceananigans
using Oceananigans.Units

using CairoMakie
using CUDA
using Printf
using Random
using SeawaterPolynomials.TEOS10: TEOS10EquationOfState

Random.seed!(1969) # for reproducible results
Random.TaskLocalRNG()

The grid

We use 128²×64 grid points with 1 m grid spacing in the horizontal and varying spacing in the vertical, with higher resolution closer to the surface. Here we use a stretching function for the vertical nodes that maintains relatively constant vertical spacing in the mixed layer, which is desirable from a numerical standpoint:

julia
Nx = Ny = 128    # number of points in each of horizontal directions
Nz = 64          # number of points in the vertical direction

Lx = Ly = 128    # (m) domain horizontal extents
Lz = 64          # (m) domain depth

refinement = 1.2 # controls spacing near surface (higher means finer spaced)
stretching = 12  # controls rate of stretching at bottom

# Normalized height ranging from 0 to 1
h(k) = (k - 1) / Nz

# Linear near-surface generator
ζ₀(k) = 1 + (h(k) - 1) / refinement

# Bottom-intensified stretching function
Σ(k) = (1 - exp(-stretching * h(k))) / (1 - exp(-stretching))

# Generating function
z_interfaces(k) = Lz * (ζ₀(k) * Σ(k) - 1)

grid = RectilinearGrid(GPU(),
                       size = (Nx, Nx, Nz),
                       x = (0, Lx),
                       y = (0, Ly),
                       z = z_interfaces)
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] variably spaced with min(Δz)=0.833413, max(Δz)=1.96618

We plot vertical spacing versus depth to inspect the prescribed grid stretching:

julia
fig = Figure(size=(1200, 800))
ax = Axis(fig[1, 1], ylabel = "z (m)", xlabel = "Vertical spacing (m)")

lines!(ax, zspacings(grid, Center()))
scatter!(ax, zspacings(grid, Center()))

fig

Buoyancy that depends on temperature and salinity

We use the SeawaterBuoyancy model with the TEOS10 equation of state,

julia
ρₒ = 1026 # kg m⁻³, average density at the surface of the world ocean
equation_of_state = TEOS10EquationOfState(reference_density=ρₒ)
buoyancy = SeawaterBuoyancy(; equation_of_state)
SeawaterBuoyancy{Float64}:
├── gravitational_acceleration: 9.80665
└── equation_of_state: BoussinesqEquationOfState{Float64}

Boundary conditions

We calculate the surface temperature flux associated with surface cooling of 200 W m⁻², reference density ρₒ, and heat capacity cᴾ,

julia
Q = 200   # W m⁻², surface _heat_ flux
cᴾ = 3991 # J K⁻¹ kg⁻¹, typical heat capacity for seawater

Jᵀ = Q / (ρₒ * cᴾ) # K m s⁻¹, surface _temperature_ flux
4.884283985946938e-5

Finally, we impose a temperature gradient dTdz both initially (see "Initial conditions" section below) and at the bottom of the domain, culminating in the boundary conditions on temperature,

julia
dTdz = 0.01 # K m⁻¹

T_bcs = FieldBoundaryConditions(top = FluxBoundaryCondition(Jᵀ),
                                bottom = GradientBoundaryCondition(dTdz))
Oceananigans.FieldBoundaryConditions, with boundary conditions
├── west: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── east: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── south: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── north: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── bottom: GradientBoundaryCondition: 0.01
├── top: FluxBoundaryCondition: 4.88428e-5
└── immersed: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)

Note that a positive temperature flux at the surface of the ocean implies cooling. This is because a positive temperature flux implies that temperature is fluxed upwards, out of the ocean.

For the velocity field, we imagine a wind blowing over the ocean surface with an average velocity at 10 meters u₁₀, and use a drag coefficient cᴰ to estimate the kinematic stress (that is, stress divided by density) exerted by the wind on the ocean:

julia
u₁₀ = 10  # m s⁻¹, average wind velocity 10 meters above the ocean
cᴰ = 2e-3 # dimensionless drag coefficient
ρₐ = 1.2  # kg m⁻³, approximate average density of air at sea-level
τx = - ρₐ / ρₒ * cᴰ * u₁₀ * abs(u₁₀) # m² s⁻²
-0.00023391812865497074

The boundary conditions on u are thus

julia
u_bcs = 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: -0.000233918
└── immersed: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)

For salinity, S, we impose an evaporative flux of the form

julia
@inline(x, y, t, S, evaporation_rate) = - evaporation_rate * S # [salinity unit] m s⁻¹

where S is salinity. We use an evaporation rate of 1 millimeter per hour,

julia
evaporation_rate = 1e-3 / hour # m s⁻¹
2.7777777777777776e-7

We build the Flux evaporation BoundaryCondition with the function , indicating that depends on salinity S and passing the parameter evaporation_rate,

julia
evaporation_bc = FluxBoundaryCondition(Jˢ, field_dependencies=:S, parameters=evaporation_rate)
FluxBoundaryCondition: ContinuousBoundaryFunction Jˢ at (Nothing, Nothing, Nothing)

The full salinity boundary conditions are

julia
S_bcs = FieldBoundaryConditions(top=evaporation_bc)
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: ContinuousBoundaryFunction Jˢ at (Nothing, Nothing, Nothing)
└── immersed: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)

Model instantiation

We fill in the final details of the model here, i.e., Coriolis forces, advection scheme, and use the DynamicSmagorinsky closure for large eddy simulation. The effect of both the WENO advection scheme and DynamicSmagorinsky turbulence closure is to dissipate variance at the grid scale. In the context of large eddy simulation, this dissipation may be interpreted as approximating a forward cascade of kinetic energy from resolved motions into motions that are smaller than the grid scale and not explicitly resolved. Note that dissipation of grid-scale variance can be achieved by the WENO advection scheme alone with closure = nothing. Typically, using an explicit closure = DynamicSmagorinsky() will produce stronger dissipation of kinetic energy; whether or not this leads to a higher quality numerical solution depends on the context. An explicit closure = DynamicSmagorinsky() is also useful for diagnosing the kinetic energy dissipation rate (the dissipation rate associated with WENO advection can be computed in principle, but is challenging and relatively computationally intensive).

In the DynamicSmagorinsky closure used below, a dynamic, multi-scale method is used to estimate the Smagorinsky coefficient at every point in time and space. Specifically, we use an algorithm that assumes the coefficient does not depend on the spatial scale of the implicit filter that separates the "true", underlying, and unresolved flow from the "filtered", or computed flow. This implementation corresponds to the "scale-invariant" formulation described by Bou-Zeid et al. (2005).

julia
model = NonhydrostaticModel(grid; buoyancy,
                            advection = WENO(order=7),
                            tracers = (:T, :S),
                            coriolis = FPlane(f=1e-4),
                            closure = DynamicSmagorinsky(),
                            boundary_conditions = (u=u_bcs, T=T_bcs, S=S_bcs))
NonhydrostaticModel{CUDAGPU, RectilinearGrid}(time = 0 seconds, iteration = 0)
├── grid: 128×128×64 RectilinearGrid{Float64, Periodic, Periodic, Bounded} on CUDAGPU with 4×4×4 halo
├── timestepper: RungeKutta3TimeStepper
├── advection scheme: WENO{4, Float64, Float32}(order=7)
├── tracers: (T, S)
├── closure: Smagorinsky with coefficient = DynamicCoefficient(averaging = LagrangianAveraging(), schedule = IterationInterval(1, 0)), Pr=(T = 1.0, S = 1.0)
├── buoyancy: SeawaterBuoyancy with g=9.80665 and BoussinesqEquationOfState{Float64} with ĝ = NegativeZDirection()
└── coriolis: FPlane{Float64}(f=0.0001)

Note: To use the (constant-coefficient) Smagorinsky-Lilly turbulence closure rather than DynamicSmagorinsky, use closure = SmagorinskyLilly() in the model constructor.

Initial conditions

Our initial condition for temperature consists of a linear stratification superposed with random noise damped at the walls, while our initial condition for velocity consists only of random noise.

julia
# Random noise damped at top and bottom
Ξ(z) = randn() * z / model.grid.Lz * (1 + z / model.grid.Lz) # noise

# Temperature initial condition: a stable density gradient with random noise superposed.
Tᵢ(x, y, z) = 20 + dTdz * z + dTdz * model.grid.Lz * 2e-6 * Ξ(z)

# Velocity initial condition: random noise scaled by the friction velocity.
uᵢ(x, y, z) = sqrt(abs(τx)) * 1e-3 * Ξ(z)

# `set!` the `model` fields using functions or constants:
set!(model, u=uᵢ, w=uᵢ, T=Tᵢ, S=35)

Setting up a simulation

We set-up a simulation with an initial time-step of 10 seconds that stops at 2 hours, with adaptive time-stepping and progress printing.

julia
simulation = Simulation(model, Δt=10, stop_time=2hours)
Simulation of NonhydrostaticModel{CUDAGPU, RectilinearGrid}(time = 0 seconds, iteration = 0)
├── Next time step: 10 seconds
├── run_wall_time: 0 seconds
├── run_wall_time / iteration: NaN days
├── stop_time: 2 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

The TimeStepWizard helps ensure stable time-stepping with a Courant-Freidrichs-Lewy (CFL) number of 0.7.

julia
conjure_time_step_wizard!(simulation, cfl=0.7)

Nice progress messaging is helpful:

julia
# Print a progress message
progress_message(sim) = @printf("Iteration: %04d, time: %s, Δt: %s, max(|w|) = %.1e ms⁻¹, wall time: %s\n",
                                iteration(sim), prettytime(sim), prettytime(sim.Δt),
                                maximum(abs, sim.model.velocities.w), prettytime(sim.run_wall_time))

add_callback!(simulation, progress_message, IterationInterval(200))

We then set up the simulation:

Output

We use the JLD2Writer to save slices of the velocity fields, tracer fields, and eddy diffusivities. The prefix keyword argument to JLD2Writer indicates that output will be saved in ocean_wind_mixing_and_convection.jld2.

julia
# Create a NamedTuple with eddy viscosity
eddy_viscosity = (; νₑ = model.closure_fields.νₑ)

filename = "ocean_wind_mixing_and_convection"

simulation.output_writers[:slices] =
    JLD2Writer(model, merge(model.velocities, model.tracers, eddy_viscosity),
               filename = filename * ".jld2",
               indices = (:, grid.Ny/2, :),
               schedule = TimeInterval(1minute),
               overwrite_existing = true)
JLD2Writer scheduled on TimeInterval(1 minute):
├── filepath: ocean_wind_mixing_and_convection.jld2
├── 6 outputs: (u, v, w, T, S, νₑ)
├── array_type: Array{Float32}
├── including: [:grid, :coriolis, :buoyancy, :closure]
├── file_splitting: NoFileSplitting
└── file size: 0 bytes (file not yet created)

We're ready:

julia
run!(simulation)
[ Info: Initializing simulation...
Iteration: 0000, time: 0 seconds, Δt: 11 seconds, max(|w|) = 1.4e-05 ms⁻¹, wall time: 0 seconds
[ Info:     ... simulation initialization complete (14.252 seconds)
[ Info: Executing initial time step...
[ Info:     ... initial time step complete (6.416 seconds).
Iteration: 0200, time: 15.999 minutes, Δt: 2.483 seconds, max(|w|) = 8.2e-06 ms⁻¹, wall time: 40.136 seconds
Iteration: 0400, time: 22.703 minutes, Δt: 1.713 seconds, max(|w|) = 1.4e-04 ms⁻¹, wall time: 1.006 minutes
Iteration: 0600, time: 27.541 minutes, Δt: 1.132 seconds, max(|w|) = 1.3e-02 ms⁻¹, wall time: 1.358 minutes
Iteration: 0800, time: 30.453 minutes, Δt: 780.078 ms, max(|w|) = 1.8e-01 ms⁻¹, wall time: 1.650 minutes
Iteration: 1000, time: 33.570 minutes, Δt: 1.312 seconds, max(|w|) = 1.5e-01 ms⁻¹, wall time: 1.946 minutes
Iteration: 1200, time: 38.962 minutes, Δt: 1.842 seconds, max(|w|) = 8.9e-02 ms⁻¹, wall time: 2.241 minutes
Iteration: 1400, time: 45.406 minutes, Δt: 1.886 seconds, max(|w|) = 7.5e-02 ms⁻¹, wall time: 2.563 minutes
Iteration: 1600, time: 51.695 minutes, Δt: 1.815 seconds, max(|w|) = 6.6e-02 ms⁻¹, wall time: 2.886 minutes
Iteration: 1800, time: 57.800 minutes, Δt: 1.927 seconds, max(|w|) = 7.0e-02 ms⁻¹, wall time: 3.219 minutes
Iteration: 2000, time: 1.060 hours, Δt: 1.812 seconds, max(|w|) = 7.8e-02 ms⁻¹, wall time: 3.576 minutes
Iteration: 2200, time: 1.151 hours, Δt: 1.703 seconds, max(|w|) = 7.6e-02 ms⁻¹, wall time: 3.868 minutes
Iteration: 2400, time: 1.246 hours, Δt: 1.615 seconds, max(|w|) = 6.6e-02 ms⁻¹, wall time: 4.156 minutes
Iteration: 2600, time: 1.335 hours, Δt: 1.583 seconds, max(|w|) = 6.9e-02 ms⁻¹, wall time: 4.459 minutes
Iteration: 2800, time: 1.425 hours, Δt: 1.699 seconds, max(|w|) = 7.1e-02 ms⁻¹, wall time: 4.742 minutes
Iteration: 3000, time: 1.516 hours, Δt: 1.660 seconds, max(|w|) = 7.3e-02 ms⁻¹, wall time: 5.024 minutes
Iteration: 3200, time: 1.604 hours, Δt: 1.579 seconds, max(|w|) = 7.5e-02 ms⁻¹, wall time: 5.332 minutes
Iteration: 3400, time: 1.691 hours, Δt: 1.586 seconds, max(|w|) = 7.3e-02 ms⁻¹, wall time: 5.643 minutes
Iteration: 3600, time: 1.776 hours, Δt: 1.564 seconds, max(|w|) = 7.8e-02 ms⁻¹, wall time: 5.997 minutes
Iteration: 3800, time: 1.860 hours, Δt: 1.438 seconds, max(|w|) = 7.7e-02 ms⁻¹, wall time: 6.346 minutes
Iteration: 4000, time: 1.940 hours, Δt: 1.511 seconds, max(|w|) = 8.0e-02 ms⁻¹, wall time: 6.730 minutes
[ Info: Simulation is stopping after running for 6.935 minutes.
[ Info: Simulation time 2 hours equals or exceeds stop time 2 hours.

Turbulence visualization

We animate the data saved in ocean_wind_mixing_and_convection.jld2. We prepare for animating the flow by loading the data into FieldTimeSeries and defining functions for computing colorbar limits.

julia
filepath = filename * ".jld2"

time_series = (w = FieldTimeSeries(filepath, "w"),
               T = FieldTimeSeries(filepath, "T"),
               S = FieldTimeSeries(filepath, "S"),
               νₑ = FieldTimeSeries(filepath, "νₑ"))
(w = 128×1×65×121 FieldTimeSeries{InMemory} located at (Center, Center, Face) of w at ocean_wind_mixing_and_convection.jld2
├── grid: 128×128×64 RectilinearGrid{Float64, Periodic, Periodic, Bounded} on CPU with 4×4×4 halo
├── indices: (:, 64:64, :)
├── time_indexing: Linear()
├── backend: InMemory()
├── path: ocean_wind_mixing_and_convection.jld2
├── name: w
└── data: 136×1×73×121 OffsetArray(::Array{Float64, 4}, -3:132, 64:64, -3:69, 1:121) with eltype Float64 with indices -3:132×64:64×-3:69×1:121
    └── max=0.0676256, min=-0.0969892, mean=-9.62242e-6, T = 128×1×64×121 FieldTimeSeries{InMemory} located at (Center, Center, Center) of T at ocean_wind_mixing_and_convection.jld2
├── grid: 128×128×64 RectilinearGrid{Float64, Periodic, Periodic, Bounded} on CPU with 4×4×4 halo
├── indices: (:, 64:64, :)
├── time_indexing: Linear()
├── backend: InMemory()
├── path: ocean_wind_mixing_and_convection.jld2
├── name: T
└── data: 136×1×72×121 OffsetArray(::Array{Float64, 4}, -3:132, 64:64, -3:68, 1:121) with eltype Float64 with indices -3:132×64:64×-3:68×1:121
    └── max=19.9983, min=0.0, mean=18.072, S = 128×1×64×121 FieldTimeSeries{InMemory} located at (Center, Center, Center) of S at ocean_wind_mixing_and_convection.jld2
├── grid: 128×128×64 RectilinearGrid{Float64, Periodic, Periodic, Bounded} on CPU with 4×4×4 halo
├── indices: (:, 64:64, :)
├── time_indexing: Linear()
├── backend: InMemory()
├── path: ocean_wind_mixing_and_convection.jld2
├── name: S
└── data: 136×1×72×121 OffsetArray(::Array{Float64, 4}, -3:132, 64:64, -3:68, 1:121) with eltype Float64 with indices -3:132×64:64×-3:68×1:121
    └── max=35.0293, min=0.0, mean=32.084, νₑ = 128×1×64×121 FieldTimeSeries{InMemory} located at (Center, Center, Center) of νₑ at ocean_wind_mixing_and_convection.jld2
├── grid: 128×128×64 RectilinearGrid{Float64, Periodic, Periodic, Bounded} on CPU with 4×4×4 halo
├── indices: (:, 64:64, :)
├── time_indexing: Linear()
├── backend: InMemory()
├── path: ocean_wind_mixing_and_convection.jld2
├── name: νₑ
└── data: 136×1×72×121 OffsetArray(::Array{Float64, 4}, -3:132, 64:64, -3:68, 1:121) with eltype Float64 with indices -3:132×64:64×-3:68×1:121
    └── max=0.0106723, min=0.0, mean=5.3387e-5)

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.

julia
times = time_series.w.times
n = Observable(length(times))

 wₙ = @lift time_series.w[$n]
 Tₙ = @lift time_series.T[$n]
 Sₙ = @lift time_series.S[$n]
νₑₙ = @lift time_series.νₑ[$n]

fig = Figure(size = (1800, 900))

axis_kwargs = (xlabel="x (m)",
               ylabel="z (m)",
               aspect = AxisAspect(grid.Lx/grid.Lz),
               limits = ((0, grid.Lx), (-grid.Lz, 0)))

ax_w  = Axis(fig[2, 1]; title = "Vertical velocity", axis_kwargs...)
ax_T  = Axis(fig[2, 3]; title = "Temperature", axis_kwargs...)
ax_S  = Axis(fig[3, 1]; title = "Salinity", axis_kwargs...)
ax_νₑ = Axis(fig[3, 3]; title = "Eddy viscocity", axis_kwargs...)

title = @lift @sprintf("t = %s", prettytime(times[$n]))

 wlims = (-0.05, 0.05)
 Tlims = (19.7, 19.99)
 Slims = (35, 35.005)
νₑlims = (1e-6, 5e-3)

hm_w = heatmap!(ax_w, wₙ; colormap = :balance, colorrange = wlims)
Colorbar(fig[2, 2], hm_w; label = "m s⁻¹")

hm_T = heatmap!(ax_T, Tₙ; colormap = :thermal, colorrange = Tlims)
Colorbar(fig[2, 4], hm_T; label = "ᵒC")

hm_S = heatmap!(ax_S, Sₙ; colormap = :haline, colorrange = Slims)
Colorbar(fig[3, 2], hm_S; label = "g / kg")

hm_νₑ = heatmap!(ax_νₑ, νₑₙ; colormap = :thermal, colorrange = νₑlims)
Colorbar(fig[3, 4], hm_νₑ; label = "m s⁻²")

fig[1, 1:4] = Label(fig, title, fontsize=24, tellwidth=false)

fig

And now record a movie. We start the animation at   minutes since things are pretty boring till then:

julia
intro = searchsortedfirst(times, 10minutes)
frames = intro:length(times)

@info "Making a motion picture of ocean wind mixing and convection..."

CairoMakie.record(fig, filename * ".mp4", frames, framerate=8) do i
    n[] = i
end
[ Info: Making a motion picture of ocean wind mixing and convection...


Julia version and environment information

This example was executed with the following version of Julia:

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
  JULIA_PROJECT = /var/lib/buildkite-agent/Oceananigans.jl-30072/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-30072/docs/.CondaPkg/.pixi/envs/default/bin/python
  JULIA_DEBUG = Literate

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

julia
import Pkg
Pkg.status()
Status `~/Oceananigans.jl-30072/docs/Project.toml`
  [79e6a3ab] Adapt v4.5.0
  [052768ef] CUDA v5.11.0
  [13f3f980] CairoMakie v0.15.9
  [e30172f5] Documenter v1.17.0
  [daee34ce] DocumenterCitations v1.4.1
  [4710194d] DocumenterVitepress v0.3.2
  [033835bb] JLD2 v0.6.3
  [63c18a36] KernelAbstractions v0.9.40
  [98b081ad] Literate v2.21.0
  [da04e1cc] MPI v0.20.24
  [85f8d34a] NCDatasets v0.14.13
  [9e8cae18] Oceananigans v0.105.5 `..`
  [f27b6e38] Polynomials v4.1.1
  [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.