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: 15.079 seconds, max(u): (0.000e+00, 0.000e+00, 0.000e+00) m/s, next Δt: 20 minutes
[ Info: ... simulation initialization complete (15.593 seconds)
[ Info: Executing initial time step...
[ Info: ... initial time step complete (4.842 seconds).
[06.94%] i: 100, t: 1.389 days, wall time: 17.659 seconds, max(u): (1.229e-01, 1.371e-01, 1.728e-03) m/s, next Δt: 20 minutes
[13.89%] i: 200, t: 2.778 days, wall time: 1.071 seconds, max(u): (2.176e-01, 1.931e-01, 2.036e-03) m/s, next Δt: 20 minutes
[20.83%] i: 300, t: 4.167 days, wall time: 939.047 ms, max(u): (2.826e-01, 2.563e-01, 1.964e-03) m/s, next Δt: 20 minutes
[27.78%] i: 400, t: 5.556 days, wall time: 903.412 ms, max(u): (3.672e-01, 3.637e-01, 1.912e-03) m/s, next Δt: 20 minutes
[34.72%] i: 500, t: 6.944 days, wall time: 1.075 seconds, max(u): (4.608e-01, 5.444e-01, 2.224e-03) m/s, next Δt: 20 minutes
[41.67%] i: 600, t: 8.333 days, wall time: 1.171 seconds, max(u): (5.839e-01, 8.443e-01, 3.001e-03) m/s, next Δt: 20 minutes
[48.61%] i: 700, t: 9.722 days, wall time: 975.365 ms, max(u): (9.384e-01, 1.204e+00, 4.191e-03) m/s, next Δt: 20 minutes
[55.56%] i: 800, t: 11.111 days, wall time: 1.096 seconds, max(u): (1.202e+00, 1.203e+00, 4.917e-03) m/s, next Δt: 20 minutes
[62.50%] i: 900, t: 12.500 days, wall time: 929.074 ms, max(u): (1.232e+00, 1.232e+00, 4.738e-03) m/s, next Δt: 20 minutes
[69.44%] i: 1000, t: 13.889 days, wall time: 978.060 ms, max(u): (1.366e+00, 1.186e+00, 4.960e-03) m/s, next Δt: 20 minutes
[76.39%] i: 1100, t: 15.278 days, wall time: 952.698 ms, max(u): (1.228e+00, 1.040e+00, 4.484e-03) m/s, next Δt: 20 minutes
[83.33%] i: 1200, t: 16.667 days, wall time: 949.540 ms, max(u): (1.341e+00, 1.253e+00, 3.361e-03) m/s, next Δt: 20 minutes
[90.28%] i: 1300, t: 18.056 days, wall time: 986.221 ms, max(u): (1.256e+00, 1.372e+00, 3.279e-03) m/s, next Δt: 20 minutes
[97.22%] i: 1400, t: 19.444 days, wall time: 962.916 ms, max(u): (1.387e+00, 1.214e+00, 3.492e-03) m/s, next Δt: 20 minutes
[ Info: Simulation is stopping after running for 34.972 seconds.
[ Info: Simulation time 20 days equals or exceeds stop time 20 days.
[ Info: Simulation completed in 35.013 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.009366166777908802 -0.008113288320600986 -0.006850333884358406 -0.005609706975519657 -0.004384476225823164 -0.003114648861810565 -0.0018515287665650249 -0.0006155974115245044; -0.009384664706885815 -0.008122437633574009 -0.006887455005198717 -0.005628624930977821 -0.004375547636300325 -0.00312762800604105 -0.0018791548209264874 -0.0006433126982301474; -0.009396584704518318 -0.008110160939395428 -0.006881497800350189 -0.0056485445238649845 -0.0043703289702534676 -0.00312466430477798 -0.0018486323533579707 -0.0006149145192466676; -0.009391014464199543 -0.00811613816767931 -0.006878576707094908 -0.005611521657556295 -0.004390213172882795 -0.003152332967147231 -0.0018383272690698504 -0.0006129884859547019; -0.009341083467006683 -0.008078400045633316 -0.006871701218187809 -0.005620880052447319 -0.0043721385300159454 -0.0031326233875006437 -0.001890357118099928 -0.0006401976570487022; -0.009384706616401672 -0.008135131560266018 -0.006871162448078394 -0.005590117536485195 -0.004365395754575729 -0.003112931502982974 -0.0018645607633516192 -0.0006397790275514126; -0.009361064992845058 -0.00811931025236845 -0.006876045372337103 -0.005615768022835255 -0.00433762650936842 -0.0031172637827694416 -0.0018804306164383888 -0.0006333930068649352; -0.009369824081659317 -0.00810187216848135 -0.006868875119835138 -0.00561480550095439 -0.004374428652226925 -0.0031134705059230328 -0.0018781623803079128 -0.0006153315771371126; -0.009368610568344593 -0.008121367543935776 -0.006874233949929476 -0.00561577919870615 -0.004383185412734747 -0.0031417650170624256 -0.0018698893254622817 -0.0006328733870759606; -0.009373823180794716 -0.008111921139061451 -0.006886753719300032 -0.005629670340567827 -0.0043888031505048275 -0.00311722862534225 -0.0018909528153017163 -0.000607053458224982; -0.009353553876280785 -0.008103744126856327 -0.006889615207910538 -0.005622982047498226 -0.004369155038148165 -0.003111846512183547 -0.0018796109361574054 -0.0006410802598111331; -0.009359938092529774 -0.008147044107317924 -0.006875252351164818 -0.005629188846796751 -0.004369044676423073 -0.0031445336062461138 -0.0018841794226318598 -0.0006471001543104649; -0.009384891018271446 -0.008124087937176228 -0.006900106091052294 -0.005640741437673569 -0.004394338931888342 -0.003131415694952011 -0.0018894159002229571 -0.0006085348431952298; -0.009382163174450397 -0.008110662922263145 -0.006866175681352615 -0.005627476144582033 -0.004389412701129913 -0.0031290212646126747 -0.0018712027231231332 -0.0005977482069283724; -0.009387758560478687 -0.008121524006128311 -0.006894206628203392 -0.005632528569549322 -0.0043566287495195866 -0.00314725236967206 -0.0018900176510214806 -0.0006170241977088153; -0.00935764517635107 -0.008118693716824055 -0.006901347544044256 -0.005624119192361832 -0.004369049333035946 -0.003111650235950947 -0.001859638374298811 -0.0006206025718711317; -0.009357430040836334 -0.008129974827170372 -0.006838330067694187 -0.005635271314531565 -0.004346258472651243 -0.0031272433698177338 -0.0018588956445455551 -0.0006509637460112572; -0.009402773343026638 -0.008121510036289692 -0.006871406454592943 -0.005629350896924734 -0.004385687876492739 -0.0031208209693431854 -0.0018824124708771706 -0.0006354177021421492; -0.009392349049448967 -0.008135776035487652 -0.006884096190333366 -0.005609624553471804 -0.004382968880236149 -0.003105826210230589 -0.0018547291401773691 -0.0006340814288705587; -0.009390417486429214 -0.008127408102154732 -0.0068680658005177975 -0.005610285326838493 -0.004365189000964165 -0.0031146591063588858 -0.0018902289448305964 -0.0006112669943831861; -0.00936703197658062 -0.008131874725222588 -0.006886731367558241 -0.00561171630397439 -0.004359742160886526 -0.0031050704419612885 -0.001877308706752956 -0.0006341421976685524; -0.00937140267342329 -0.008100437000393867 -0.006866796873509884 -0.005628880579024553 -0.004373322241008282 -0.0031349833589047194 -0.001860482501797378 -0.0006228741258382797; -0.007515241391956806 -0.006241559982299805 -0.004992206115275621 -0.003744879039004445 -0.002489120699465275 -0.0012375748483464122 1.580213756824378e-5 0.0012435894459486008; -0.005411255639046431 -0.004159173928201199 -0.0029198837000876665 -0.0016774205723777413 -0.00042378681246191263 0.0008664191118441522 0.0020634601823985577 0.003299437928944826; -0.00333658535964787 -0.002106250263750553 -0.0008285248186439276 0.0004312139062676579 0.0016546720871701837 0.00290722050704062 0.004149190615862608 0.005408566910773516; -0.0012481072917580605 1.3628269698529039e-5 0.0012601307826116681 0.002485064323991537 0.0037523191422224045 0.00498337484896183 0.0062439595349133015 0.007516916375607252; 0.0006206960533745587 0.0018627509707584977 0.0031248137820512056 0.004387972876429558 0.005639065057039261 0.006873659323900938 0.008099925704300404 0.009365375153720379; 0.0006505097262561321 0.0018630159320309758 0.003128741169348359 0.00436948798596859 0.005634090397506952 0.0068982853554189205 0.008098351769149303 0.009383220225572586; 0.0006303864647634327 0.0018910502549260855 0.0030985944904386997 0.004349756520241499 0.00560276722535491 0.006875625345855951 0.00812501274049282 0.009378370828926563; 0.0006100856116972864 0.001854105619713664 0.003143982496112585 0.004356110002845526 0.005608057137578726 0.006874545477330685 0.008146940730512142 0.00938379019498825; 0.0006400135462172329 0.0018578843446448445 0.0031261108815670013 0.00436933571472764 0.005621237680315971 0.006892910692840815 0.008134467527270317 0.009385867975652218; 0.0006608911789953709 0.001876653404906392 0.0031131221912801266 0.004379985388368368 0.005643161945044994 0.006879163905978203 0.008133341558277607 0.009374232962727547; 0.0006183107616379857 0.001881520845927298 0.00311883888207376 0.0043779052793979645 0.005607070401310921 0.006841429974883795 0.008097752928733826 0.009385332465171814; 0.0006333138444460928 0.0018676903564482927 0.0031209576409310102 0.004397485870867968 0.005619334056973457 0.006872375961393118 0.00814428273588419 0.009391029365360737; 0.0006065453635528684 0.0018647651886567473 0.0031339661218225956 0.004360990133136511 0.005622191354632378 0.00685889134183526 0.008113620802760124 0.009388923645019531; 0.0006035068654455245 0.0018704115645959973 0.003130404045805335 0.004366676788777113 0.005599881522357464 0.006883484311401844 0.00814264826476574 0.009398485533893108; 0.0006541155162267387 0.0018846872262656689 0.0031582193914800882 0.004357491619884968 0.005617012735456228 0.006887077819555998 0.00811850093305111 0.009365418925881386; 0.0006463109166361392 0.0018657464534044266 0.00311297201551497 0.004349946975708008 0.005618355702608824 0.006882279645651579 0.008106335997581482 0.009363367222249508; 0.0006095965509302914 0.001859858981333673 0.0031256566289812326 0.004368705675005913 0.005611418280750513 0.006868628319352865 0.008102534338831902 0.009376336820423603; 0.0006302096298895776 0.001895390567369759 0.003123462200164795 0.004383897874504328 0.00562366284430027 0.006854724604636431 0.008114197291433811 0.0093615110963583; 0.0006290043238550425 0.0018695819890126586 0.0031090988777577877 0.0043684388510882854 0.005655643530189991 0.006894041318446398 0.008121988736093044 0.00936682429164648; 0.0006177318282425404 0.0018802762497216463 0.003107448574155569 0.004367568995803595 0.005626602564007044 0.006885837763547897 0.008116130717098713 0.00937708467245102; 0.0005986134056001902 0.0018600597977638245 0.0031378401909023523 0.004378729499876499 0.005615090485662222 0.00689623923972249 0.008115118369460106 0.009364108555018902; 0.00060264952480793 0.0018920168513432145 0.0031283858697861433 0.004382545594125986 0.005634980276226997 0.006876807659864426 0.0081421984359622 0.009349229745566845; 0.000652784772682935 0.0018816505325958133 0.0031180756632238626 0.004360899329185486 0.0056594242341816425 0.006885460112243891 0.008102632127702236 0.00935768336057663; 0.0006142375059425831 0.0018894337117671967 0.003124625189229846 0.004411650355905294 0.005600220989435911 0.006861357484012842 0.00810643658041954 0.009397784247994423; 0.0006279056542553008 0.001855905749835074 0.003136591287329793 0.004396364558488131 0.005627815145999193 0.006886558141559362 0.00814265850931406 0.009362281300127506; 0.0006343889399431646 0.0018635543528944254 0.0031536766327917576 0.004397151060402393 0.0056061106733977795 0.006888328120112419 0.008131349459290504 0.00934994500130415])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-29683/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-29683/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-29683/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.1 `..`
[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.