Baroclinic adjustment
In this example, we simulate the evolution and equilibration of a baroclinically unstable front.
Install dependencies
First let's make sure we have all required packages installed.
using Pkg
pkg"add Oceananigans, CairoMakie"using Oceananigans
using Oceananigans.UnitsGrid
We use a three-dimensional channel that is periodic in the x direction:
Lx = 1000kilometers # east-west extent [m]
Ly = 1000kilometers # north-south extent [m]
Lz = 1kilometers # depth [m]
grid = RectilinearGrid(size = (48, 48, 8),
x = (0, Lx),
y = (-Ly/2, Ly/2),
z = (-Lz, 0),
topology = (Periodic, Bounded, Bounded))48×48×8 RectilinearGrid{Float64, Periodic, Bounded, Bounded} on CPU with 3×3×3 halo
├── Periodic x ∈ [0.0, 1.0e6) regularly spaced with Δx=20833.3
├── Bounded y ∈ [-500000.0, 500000.0] regularly spaced with Δy=20833.3
└── Bounded z ∈ [-1000.0, 0.0] regularly spaced with Δz=125.0Model
We built a HydrostaticFreeSurfaceModel with an ImplicitFreeSurface solver. Regarding Coriolis, we use a beta-plane centered at 45° South.
model = HydrostaticFreeSurfaceModel(grid;
coriolis = BetaPlane(latitude = -45),
buoyancy = BuoyancyTracer(),
tracers = :b,
momentum_advection = WENO(),
tracer_advection = WENO())HydrostaticFreeSurfaceModel{CPU, RectilinearGrid}(time = 0 seconds, iteration = 0)
├── grid: 48×48×8 RectilinearGrid{Float64, Periodic, Bounded, Bounded} on CPU with 3×3×3 halo
├── timestepper: QuasiAdamsBashforth2TimeStepper
├── tracers: b
├── closure: Nothing
├── buoyancy: BuoyancyTracer with ĝ = NegativeZDirection()
├── free surface: ImplicitFreeSurface with gravitational acceleration 9.80665 m s⁻²
│ └── solver: FFTImplicitFreeSurfaceSolver
├── advection scheme:
│ ├── momentum: WENO{3, Float64, Float32}(order=5)
│ └── b: WENO{3, Float64, Float32}(order=5)
├── vertical_coordinate: ZCoordinate
└── coriolis: BetaPlane{Float64}We start our simulation from rest with a baroclinically unstable buoyancy distribution. We use ramp(y, Δy), defined below, to specify a front with width Δy and horizontal buoyancy gradient M². We impose the front on top of a vertical buoyancy gradient N² and a bit of noise.
"""
ramp(y, Δy)
Linear ramp from 0 to 1 between -Δy/2 and +Δy/2.
For example:
```
y < -Δy/2 => ramp = 0
-Δy/2 < y < -Δy/2 => ramp = y / Δy
y > Δy/2 => ramp = 1
```
"""
ramp(y, Δy) = min(max(0, y/Δy + 1/2), 1)
N² = 1e-5 # [s⁻²] buoyancy frequency / stratification
M² = 1e-7 # [s⁻²] horizontal buoyancy gradient
Δy = 100kilometers # width of the region of the front
Δb = Δy * M² # buoyancy jump associated with the front
ϵb = 1e-2 * Δb # noise amplitude
bᵢ(x, y, z) = N² * z + Δb * ramp(y, Δy) + ϵb * randn()
set!(model, b=bᵢ)Let's visualize the initial buoyancy distribution.
using CairoMakie
set_theme!(Theme(fontsize = 20))
# Build coordinates with units of kilometers
x, y, z = 1e-3 .* nodes(grid, (Center(), Center(), Center()))
b = model.tracers.b
fig, ax, hm = heatmap(view(b, 1, :, :),
colormap = :deep,
axis = (xlabel = "y [km]",
ylabel = "z [km]",
title = "b(x=0, y, z, t=0)",
titlesize = 24))
Colorbar(fig[1, 2], hm, label = "[m s⁻²]")
figSimulation
Now let's build a Simulation.
simulation = Simulation(model, Δt=20minutes, stop_time=20days)Simulation of HydrostaticFreeSurfaceModel{CPU, RectilinearGrid}(time = 0 seconds, iteration = 0)
├── Next time step: 20 minutes
├── run_wall_time: 0 seconds
├── run_wall_time / iteration: NaN days
├── stop_time: 20 days
├── 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 entriesWe add a TimeStepWizard callback to adapt the simulation's time-step,
conjure_time_step_wizard!(simulation, IterationInterval(20), cfl=0.2, max_Δt=20minutes)Also, we add a callback to print a message about how the simulation is going,
using Printf
wall_clock = Ref(time_ns())
function print_progress(sim)
u, v, w = model.velocities
progress = 100 * (time(sim) / sim.stop_time)
elapsed = (time_ns() - wall_clock[]) / 1e9
@printf("[%05.2f%%] i: %d, t: %s, wall time: %s, max(u): (%6.3e, %6.3e, %6.3e) m/s, next Δt: %s\n",
progress, iteration(sim), prettytime(sim), prettytime(elapsed),
maximum(abs, u), maximum(abs, v), maximum(abs, w), prettytime(sim.Δt))
wall_clock[] = time_ns()
return nothing
end
add_callback!(simulation, print_progress, IterationInterval(100))Diagnostics/Output
Here, we save the buoyancy,
u, v, w = model.velocities
ζ = ∂x(v) - ∂y(u)
B = Average(b, dims=1)
U = Average(u, dims=1)
V = Average(v, dims=1)
filename = "baroclinic_adjustment"
save_fields_interval = 0.5day
slicers = (east = (grid.Nx, :, :),
north = (:, grid.Ny, :),
bottom = (:, :, 1),
top = (:, :, grid.Nz))
for side in keys(slicers)
indices = slicers[side]
simulation.output_writers[side] = JLD2Writer(model, (; b, ζ);
filename = filename * "_$(side)_slice",
schedule = TimeInterval(save_fields_interval),
overwrite_existing = true,
indices)
end
simulation.output_writers[:zonal] = JLD2Writer(model, (; b=B, u=U, v=V);
filename = filename * "_zonal_average",
schedule = TimeInterval(save_fields_interval),
overwrite_existing = true)JLD2Writer scheduled on TimeInterval(12 hours):
├── filepath: baroclinic_adjustment_zonal_average.jld2
├── 3 outputs: (b, u, v)
├── array_type: Array{Float32}
├── including: [:grid, :coriolis, :buoyancy, :closure]
├── file_splitting: NoFileSplitting
└── file size: 0 bytes (file not yet created)Now we're ready to run.
@info "Running the simulation..."
run!(simulation)
@info "Simulation completed in " * prettytime(simulation.run_wall_time)[ Info: Running the simulation...
[ Info: Initializing simulation...
[00.00%] i: 0, t: 0 seconds, wall time: 10.761 seconds, max(u): (0.000e+00, 0.000e+00, 0.000e+00) m/s, next Δt: 20 minutes
[ Info: ... simulation initialization complete (16.329 seconds)
[ Info: Executing initial time step...
[ Info: ... initial time step complete (6.115 seconds).
[06.94%] i: 100, t: 1.389 days, wall time: 20.868 seconds, max(u): (1.232e-01, 1.309e-01, 1.687e-03) m/s, next Δt: 20 minutes
[13.89%] i: 200, t: 2.778 days, wall time: 1.208 seconds, max(u): (2.100e-01, 1.940e-01, 1.981e-03) m/s, next Δt: 20 minutes
[20.83%] i: 300, t: 4.167 days, wall time: 1.511 seconds, max(u): (2.999e-01, 2.697e-01, 1.878e-03) m/s, next Δt: 20 minutes
[27.78%] i: 400, t: 5.556 days, wall time: 1.702 seconds, max(u): (3.775e-01, 4.012e-01, 1.909e-03) m/s, next Δt: 20 minutes
[34.72%] i: 500, t: 6.944 days, wall time: 1.601 seconds, max(u): (5.065e-01, 6.253e-01, 2.294e-03) m/s, next Δt: 20 minutes
[41.67%] i: 600, t: 8.333 days, wall time: 1.607 seconds, max(u): (6.766e-01, 1.061e+00, 3.016e-03) m/s, next Δt: 20 minutes
[48.61%] i: 700, t: 9.722 days, wall time: 1.565 seconds, max(u): (1.023e+00, 1.325e+00, 4.939e-03) m/s, next Δt: 20 minutes
[55.56%] i: 800, t: 11.111 days, wall time: 1.703 seconds, max(u): (1.388e+00, 1.173e+00, 4.415e-03) m/s, next Δt: 20 minutes
[62.50%] i: 900, t: 12.500 days, wall time: 1.405 seconds, max(u): (1.350e+00, 1.297e+00, 5.283e-03) m/s, next Δt: 20 minutes
[69.44%] i: 1000, t: 13.889 days, wall time: 1.605 seconds, max(u): (1.362e+00, 1.062e+00, 4.531e-03) m/s, next Δt: 20 minutes
[76.39%] i: 1100, t: 15.278 days, wall time: 1.534 seconds, max(u): (1.434e+00, 1.064e+00, 2.888e-03) m/s, next Δt: 20 minutes
[83.33%] i: 1200, t: 16.667 days, wall time: 1.497 seconds, max(u): (1.280e+00, 1.184e+00, 3.411e-03) m/s, next Δt: 20 minutes
[90.28%] i: 1300, t: 18.056 days, wall time: 1.903 seconds, max(u): (1.242e+00, 1.311e+00, 2.344e-03) m/s, next Δt: 20 minutes
[97.22%] i: 1400, t: 19.444 days, wall time: 1.110 seconds, max(u): (1.379e+00, 1.477e+00, 3.557e-03) m/s, next Δt: 20 minutes
[ Info: Simulation is stopping after running for 43.786 seconds.
[ Info: Simulation time 20 days equals or exceeds stop time 20 days.
[ Info: Simulation completed in 43.824 secondsVisualization
All that's left is to make a pretty movie. Actually, we make two visualizations here. First, we illustrate how to make a 3D visualization with Makie's Axis3 and Makie.surface. Then we make a movie in 2D. We use CairoMakie in this example, but note that using GLMakie is more convenient on a system with OpenGL, as figures will be displayed on the screen.
using CairoMakieThree-dimensional visualization
We load the saved buoyancy output on the top, north, and east surface as FieldTimeSerieses.
filename = "baroclinic_adjustment"
sides = keys(slicers)
slice_filenames = NamedTuple(side => filename * "_$(side)_slice.jld2" for side in sides)
b_timeserieses = (east = FieldTimeSeries(slice_filenames.east, "b"),
north = FieldTimeSeries(slice_filenames.north, "b"),
top = FieldTimeSeries(slice_filenames.top, "b"))
B_timeseries = FieldTimeSeries(filename * "_zonal_average.jld2", "b")
times = B_timeseries.times
grid = B_timeseries.grid48×48×8 RectilinearGrid{Float64, Periodic, Bounded, Bounded} on CPU with 3×3×3 halo
├── Periodic x ∈ [0.0, 1.0e6) regularly spaced with Δx=20833.3
├── Bounded y ∈ [-500000.0, 500000.0] regularly spaced with Δy=20833.3
└── Bounded z ∈ [-1000.0, 0.0] regularly spaced with Δz=125.0We build the coordinates. We rescale horizontal coordinates to kilometers.
xb, yb, zb = nodes(b_timeserieses.east)
xb = xb ./ 1e3 # convert m -> km
yb = yb ./ 1e3 # convert m -> km
Nx, Ny, Nz = size(grid)
x_xz = repeat(x, 1, Nz)
y_xz_north = y[end] * ones(Nx, Nz)
z_xz = repeat(reshape(z, 1, Nz), Nx, 1)
x_yz_east = x[end] * ones(Ny, Nz)
y_yz = repeat(y, 1, Nz)
z_yz = repeat(reshape(z, 1, Nz), grid.Ny, 1)
x_xy = x
y_xy = y
z_xy_top = z[end] * ones(grid.Nx, grid.Ny)Then we create a 3D axis. We use zonal_slice_displacement to control where the plot of the instantaneous zonal average flow is located.
fig = Figure(size = (1600, 800))
zonal_slice_displacement = 1.2
ax = Axis3(fig[2, 1],
aspect=(1, 1, 1/5),
xlabel = "x (km)",
ylabel = "y (km)",
zlabel = "z (m)",
xlabeloffset = 100,
ylabeloffset = 100,
zlabeloffset = 100,
limits = ((x[1], zonal_slice_displacement * x[end]), (y[1], y[end]), (z[1], z[end])),
elevation = 0.45,
azimuth = 6.8,
xspinesvisible = false,
zgridvisible = false,
protrusions = 40,
perspectiveness = 0.7)Axis3 with 12 plots:
┣━ Poly{Tuple{GeometryBasics.Polygon{2, Float64}}}
┣━ Poly{Tuple{GeometryBasics.Polygon{2, Float64}}}
┣━ Poly{Tuple{GeometryBasics.Polygon{2, Float64}}}
┣━ LineSegments{Tuple{Base.ReinterpretArray{Point{3, Float64}, 1, Tuple{Point{3, Float64}, Point{3, Float64}}, Vector{Tuple{Point{3, Float64}, Point{3, Float64}}}, false}}}
┣━ LineSegments{Tuple{Base.ReinterpretArray{Point{3, Float64}, 1, Tuple{Point{3, Float64}, Point{3, Float64}}, Vector{Tuple{Point{3, Float64}, Point{3, Float64}}}, false}}}
┣━ LineSegments{Tuple{Vector{Point{3, Float64}}}}
┣━ LineSegments{Tuple{Base.ReinterpretArray{Point{3, Float64}, 1, Tuple{Point{3, Float64}, Point{3, Float64}}, Vector{Tuple{Point{3, Float64}, Point{3, Float64}}}, false}}}
┣━ LineSegments{Tuple{Base.ReinterpretArray{Point{3, Float64}, 1, Tuple{Point{3, Float64}, Point{3, Float64}}, Vector{Tuple{Point{3, Float64}, Point{3, Float64}}}, false}}}
┣━ LineSegments{Tuple{Vector{Point{3, Float64}}}}
┣━ LineSegments{Tuple{Base.ReinterpretArray{Point{3, Float64}, 1, Tuple{Point{3, Float64}, Point{3, Float64}}, Vector{Tuple{Point{3, Float64}, Point{3, Float64}}}, false}}}
┣━ LineSegments{Tuple{Base.ReinterpretArray{Point{3, Float64}, 1, Tuple{Point{3, Float64}, Point{3, Float64}}, Vector{Tuple{Point{3, Float64}, Point{3, Float64}}}, false}}}
┗━ LineSegments{Tuple{Vector{Point{3, Float64}}}}We use data from the final savepoint for the 3D plot. Note that this plot can easily be animated by using Makie's Observable. To dive into Observables, check out Makie.jl's Documentation.
n = length(times)41Now let's make a 3D plot of the buoyancy and in front of it we'll use the zonally-averaged output to plot the instantaneous zonal-average of the buoyancy.
b_slices = (east = interior(b_timeserieses.east[n], 1, :, :),
north = interior(b_timeserieses.north[n], :, 1, :),
top = interior(b_timeserieses.top[n], :, :, 1))
# Zonally-averaged buoyancy
B = interior(B_timeseries[n], 1, :, :)
clims = 1.1 .* extrema(b_timeserieses.top[n][:])
kwargs = (colorrange=clims, colormap=:deep, shading=NoShading)
surface!(ax, x_yz_east, y_yz, z_yz; color = b_slices.east, kwargs...)
surface!(ax, x_xz, y_xz_north, z_xz; color = b_slices.north, kwargs...)
surface!(ax, x_xy, y_xy, z_xy_top; color = b_slices.top, kwargs...)
sf = surface!(ax, zonal_slice_displacement .* x_yz_east, y_yz, z_yz; color = B, kwargs...)
contour!(ax, y, z, B; transformation = (:yz, zonal_slice_displacement * x[end]),
levels = 15, linewidth = 2, color = :black)
Colorbar(fig[2, 2], sf, label = "m s⁻²", height = Relative(0.4), tellheight=false)
title = "Buoyancy at t = " * string(round(times[n] / day, digits=1)) * " days"
fig[1, 1:2] = Label(fig, title; fontsize = 24, tellwidth = false, padding = (0, 0, -120, 0))
rowgap!(fig.layout, 1, Relative(-0.2))
colgap!(fig.layout, 1, Relative(-0.1))
save("baroclinic_adjustment_3d.png", fig)
Two-dimensional movie
We make a 2D movie that shows buoyancy FieldTimeSeries and extract the additional coordinates we'll need for plotting
ζ_timeseries = FieldTimeSeries(slice_filenames.top, "ζ")
U_timeseries = FieldTimeSeries(filename * "_zonal_average.jld2", "u")
B_timeseries = FieldTimeSeries(filename * "_zonal_average.jld2", "b")
V_timeseries = FieldTimeSeries(filename * "_zonal_average.jld2", "v")
xζ, yζ, zζ = nodes(ζ_timeseries)
yv = ynodes(V_timeseries)
xζ = xζ ./ 1e3 # convert m -> km
yζ = yζ ./ 1e3 # convert m -> km
yv = yv ./ 1e3 # convert m -> km-500.0:20.833333333333332:500.0Next, we set up a plot with 4 panels. The top panels are large and square, while the bottom panels get a reduced aspect ratio through rowsize!.
fig = Figure(size=(1800, 1000))
axb = Axis(fig[1, 2], xlabel="x (km)", ylabel="y (km)", aspect=1)
axζ = Axis(fig[1, 3], xlabel="x (km)", ylabel="y (km)", aspect=1, yaxisposition=:right)
axu = Axis(fig[2, 2], xlabel="y (km)", ylabel="z (m)")
axv = Axis(fig[2, 3], xlabel="y (km)", ylabel="z (m)", yaxisposition=:right)
rowsize!(fig.layout, 2, Relative(0.3))To prepare a plot for animation, we index the timeseries with an Observable,
n = Observable(1)
b_top = @lift interior(b_timeserieses.top[$n], :, :, 1)
ζ_top = @lift interior(ζ_timeseries[$n], :, :, 1)
U = @lift interior(U_timeseries[$n], 1, :, :)
V = @lift interior(V_timeseries[$n], 1, :, :)
B = @lift interior(B_timeseries[$n], 1, :, :)Observable([-0.009363396093249321 -0.008139862678945065 -0.0068682231940329075 -0.0056414236314594746 -0.0043716710060834885 -0.003125683404505253 -0.0018739616498351097 -0.0006394109223037958; -0.009373249486088753 -0.008117950521409512 -0.006837823893874884 -0.0056234560906887054 -0.004367616958916187 -0.0031401291489601135 -0.0018650011625140905 -0.0006392996292561293; -0.00938014592975378 -0.008118054829537868 -0.006870056968182325 -0.00562239158898592 -0.00435091694816947 -0.003152113640680909 -0.0018768190639093518 -0.0006149793625809252; -0.009351040236651897 -0.008108476176857948 -0.006885851267725229 -0.005629741121083498 -0.0043751890771090984 -0.003129026386886835 -0.0018609624821692705 -0.000613053678534925; -0.00937472190707922 -0.008124767802655697 -0.006894745863974094 -0.005646656267344952 -0.004380083177238703 -0.003136188955977559 -0.0018770744791254401 -0.0006249095313251019; -0.009353948757052422 -0.008110033348202705 -0.006863179616630077 -0.005659361369907856 -0.004387400578707457 -0.003134204540401697 -0.0018903716700151563 -0.0006356898229569197; -0.009377053938806057 -0.00813875813037157 -0.00687563419342041 -0.005641954485327005 -0.00437787314876914 -0.0031288990285247564 -0.0018703207606449723 -0.0006322207045741379; -0.009359663352370262 -0.008117184042930603 -0.006864996161311865 -0.00561907421797514 -0.004360314458608627 -0.003138316562399268 -0.0018819704419001937 -0.0006450113723985851; -0.009381169453263283 -0.008123833686113358 -0.006870530545711517 -0.00562874972820282 -0.00437987269833684 -0.00311447586864233 -0.0018602946074679494 -0.0006255959160625935; -0.009358260780572891 -0.008110000751912594 -0.006864046677947044 -0.0056130895391106606 -0.00440028915181756 -0.003114885650575161 -0.0018785996362566948 -0.0006202419754117727; -0.009377989917993546 -0.008102466352283955 -0.006870089564472437 -0.0056180451065301895 -0.00436426093801856 -0.0031280231196433306 -0.0018700968939810991 -0.0006323032430373132; -0.009373322129249573 -0.00812856387346983 -0.006860848981887102 -0.005650026723742485 -0.004369833040982485 -0.0031367174815386534 -0.0018504568142816424 -0.0006341031403280795; -0.009353159926831722 -0.00812462717294693 -0.006896293722093105 -0.00561263132840395 -0.004357574041932821 -0.0031234899070113897 -0.001860524876974523 -0.0006410062778741121; -0.0093711381778121 -0.008117718622088432 -0.0068733468651771545 -0.005621564574539661 -0.004357350990176201 -0.003123067319393158 -0.001879263436421752 -0.0006205992540344596; -0.009366193786263466 -0.008118034340441227 -0.0068733287043869495 -0.005628818646073341 -0.004362171981483698 -0.003114443738013506 -0.001898249494843185 -0.0006375758093781769; -0.009394693188369274 -0.008111776784062386 -0.006874077022075653 -0.005627563688904047 -0.004351147450506687 -0.0031039107125252485 -0.001900918548926711 -0.0006180621567182243; -0.009400331415235996 -0.008116292767226696 -0.006887275725603104 -0.005627489183098078 -0.0043818214908242226 -0.003158321836963296 -0.001910109887830913 -0.0006138526368886232; -0.009350404143333435 -0.008124098181724548 -0.006885213777422905 -0.00562631618231535 -0.004359724465757608 -0.0031399365980178118 -0.0018612948479130864 -0.0005958896945230663; -0.00938791036605835 -0.00812780112028122 -0.006881110370159149 -0.005626122932881117 -0.004351698327809572 -0.0031198407523334026 -0.0018865464953705668 -0.0006239969516173005; -0.009364780969917774 -0.008109595626592636 -0.006871150806546211 -0.005617604590952396 -0.004422195255756378 -0.003111035330221057 -0.0018630324630066752 -0.0006123940111137927; -0.009364827536046505 -0.0081104077398777 -0.006883577443659306 -0.005641635973006487 -0.004360529128462076 -0.0031154972966760397 -0.0018842500867322087 -0.0006294181221164763; -0.00938292033970356 -0.008112561888992786 -0.006878370884805918 -0.005627170205116272 -0.0043731508776545525 -0.0031138299964368343 -0.0018560824682936072 -0.0006231318111531436; -0.007504087872803211 -0.006289151031523943 -0.004979229997843504 -0.0037350414786487818 -0.0024968197103589773 -0.0012393060605973005 -3.223632120352704e-6 0.001258022035472095; -0.0054299794137477875 -0.004165487363934517 -0.0029361415654420853 -0.001674428815022111 -0.00043063872726634145 0.0008512805798090994 0.00208842009305954 0.0033317632041871548; -0.0033307115081697702 -0.0020759545732289553 -0.0008261752664111555 0.0004277575935702771 0.001672930782660842 0.002909476635977626 0.00419971439987421 0.005440925247967243; -0.0012501721503213048 2.7590801892074523e-6 0.0012350552715361118 0.002500670962035656 0.003771738847717643 0.005032157059758902 0.006255371496081352 0.007478075101971626; 0.0005981709691695869 0.0018766614375635982 0.0031209427397698164 0.004407758824527264 0.005621056538075209 0.006870258133858442 0.008120989426970482 0.009375708177685738; 0.000639422214590013 0.001872522640042007 0.003149853553622961 0.004375042859464884 0.005597106646746397 0.006860633846372366 0.008141548372805119 0.009372951462864876; 0.0006196332396939397 0.001876325928606093 0.003114673774689436 0.004386737011373043 0.005623779725283384 0.006873657926917076 0.008110643364489079 0.00935387797653675; 0.0006085268105380237 0.001890755956992507 0.0031211518216878176 0.004366261884570122 0.005625409539788961 0.006867224350571632 0.008110337890684605 0.009380485862493515; 0.0006215738249011338 0.0018775039352476597 0.0031534619629383087 0.00437813950702548 0.005604149308055639 0.006873059086501598 0.008120720274746418 0.00937762763351202; 0.0006430236389860511 0.0018760422244668007 0.003128236858174205 0.004389587789773941 0.005623455159366131 0.0068954466842114925 0.008101468905806541 0.009379008784890175; 0.0006325285066850483 0.0018462707521393895 0.0031089368276298046 0.004381003323942423 0.005600492935627699 0.006880105938762426 0.008110850118100643 0.009406442753970623; 0.0006308112060651183 0.0018777056830003858 0.003130615223199129 0.004378307610750198 0.005612833891063929 0.006888469215482473 0.008109340444207191 0.00936870463192463; 0.0006035711849108338 0.001880659838207066 0.003120700130239129 0.004386104643344879 0.005641601514071226 0.0068795992992818356 0.008117707446217537 0.009365448728203773; 0.0006124976789578795 0.0018913710955530405 0.0031059677712619305 0.004362528678029776 0.005658277310431004 0.0068740444257855415 0.008130410686135292 0.009373783133924007; 0.000603552907705307 0.0018362632254138589 0.0031041118782013655 0.004365057218819857 0.0056227208115160465 0.006867375690490007 0.008099832572042942 0.009380631148815155; 0.0006328757735900581 0.0018778274534270167 0.003136072074994445 0.004389068111777306 0.005632603541016579 0.006854680832475424 0.008097237907350063 0.009356435388326645; 0.0006474804831668735 0.0018775522476062179 0.0031427552457898855 0.004380691330879927 0.005620751529932022 0.006896196864545345 0.008143884129822254 0.009395974688231945; 0.0005996887921355665 0.0018640629714354873 0.0031446439679712057 0.004374560434371233 0.005641044117510319 0.006863969843834639 0.008122418075799942 0.009389173239469528; 0.0006378756952472031 0.0018664579838514328 0.0030881823040544987 0.0043735746294260025 0.005621620453894138 0.006864946335554123 0.008132299408316612 0.00936818402260542; 0.0006465958431363106 0.0018440117128193378 0.0031057638116180897 0.004364512860774994 0.005619965493679047 0.006870206445455551 0.00810993555933237 0.009400051087141037; 0.0006389334448613226 0.0018558503361418843 0.003139347303658724 0.004388031084090471 0.0056155286729335785 0.006864850874990225 0.008133248426020145 0.00938528124243021; 0.0006385795422829688 0.0018740467494353652 0.0031080448534339666 0.004366607870906591 0.005613413173705339 0.006890412420034409 0.008125760592520237 0.00938008539378643; 0.0006190243293531239 0.001872627530246973 0.0031261956319212914 0.004360951948910952 0.0056329709477722645 0.006862349342554808 0.008107345551252365 0.00934386532753706; 0.0006090085371397436 0.001870727981440723 0.0031305397860705853 0.004384275525808334 0.005619593430310488 0.006863879039883614 0.0081039909273386 0.009385923855006695; 0.0006264548283070326 0.0018563024932518601 0.0031300203409045935 0.00437692878767848 0.005586809013038874 0.006867923308163881 0.00812326930463314 0.009396737441420555; 0.0006077754660509527 0.0018583819037303329 0.0031242503318935633 0.004384152591228485 0.0056193978525698185 0.006897030398249626 0.008141706697642803 0.009374754503369331])and then build our plot:
hm = heatmap!(axb, xb, yb, b_top, colorrange=(0, Δb), colormap=:thermal)
Colorbar(fig[1, 1], hm, flipaxis=false, label="Surface b(x, y) (m s⁻²)")
hm = heatmap!(axζ, xζ, yζ, ζ_top, colorrange=(-5e-5, 5e-5), colormap=:balance)
Colorbar(fig[1, 4], hm, label="Surface ζ(x, y) (s⁻¹)")
hm = heatmap!(axu, yb, zb, U; colorrange=(-5e-1, 5e-1), colormap=:balance)
Colorbar(fig[2, 1], hm, flipaxis=false, label="Zonally-averaged U(y, z) (m s⁻¹)")
contour!(axu, yb, zb, B; levels=15, color=:black)
hm = heatmap!(axv, yv, zb, V; colorrange=(-1e-1, 1e-1), colormap=:balance)
Colorbar(fig[2, 4], hm, label="Zonally-averaged V(y, z) (m s⁻¹)")
contour!(axv, yb, zb, B; levels=15, color=:black)Finally, we're ready to record the movie.
frames = 1:length(times)
record(fig, filename * ".mp4", frames, framerate=8) do i
n[] = i
endJulia 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-29586/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-29586/docs/.CondaPkg/.pixi/envs/default/bin/python
JULIA_DEBUG = LiterateThese were the top-level packages installed in the environment:
import Pkg
Pkg.status()Status `~/Oceananigans.jl-29586/docs/Project.toml`
[79e6a3ab] Adapt v4.4.0
[052768ef] CUDA v5.9.6
[13f3f980] CairoMakie v0.15.8
[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.23
[85f8d34a] NCDatasets v0.14.11
[9e8cae18] Oceananigans v0.105.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.1This page was generated using Literate.jl.