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, $b$, at the edges of our domain as well as the zonal ($x$) average of 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: 9.045 seconds, max(u): (0.000e+00, 0.000e+00, 0.000e+00) m/s, next Δt: 20 minutes
[ Info: ... simulation initialization complete (17.842 seconds)
[ Info: Executing initial time step...
[ Info: ... initial time step complete (6.063 seconds).
[06.94%] i: 100, t: 1.389 days, wall time: 21.991 seconds, max(u): (1.291e-01, 1.220e-01, 1.874e-03) m/s, next Δt: 20 minutes
[13.89%] i: 200, t: 2.778 days, wall time: 1.170 seconds, max(u): (2.109e-01, 1.963e-01, 1.957e-03) m/s, next Δt: 20 minutes
[20.83%] i: 300, t: 4.167 days, wall time: 862.718 ms, max(u): (2.930e-01, 2.657e-01, 1.858e-03) m/s, next Δt: 20 minutes
[27.78%] i: 400, t: 5.556 days, wall time: 1.005 seconds, max(u): (3.884e-01, 4.237e-01, 1.716e-03) m/s, next Δt: 20 minutes
[34.72%] i: 500, t: 6.944 days, wall time: 875.798 ms, max(u): (4.865e-01, 5.744e-01, 1.980e-03) m/s, next Δt: 20 minutes
[41.67%] i: 600, t: 8.333 days, wall time: 959.338 ms, max(u): (6.058e-01, 8.264e-01, 2.628e-03) m/s, next Δt: 20 minutes
[48.61%] i: 700, t: 9.722 days, wall time: 1.043 seconds, max(u): (8.306e-01, 1.213e+00, 3.391e-03) m/s, next Δt: 20 minutes
[55.56%] i: 800, t: 11.111 days, wall time: 1.053 seconds, max(u): (1.181e+00, 1.342e+00, 4.561e-03) m/s, next Δt: 20 minutes
[62.50%] i: 900, t: 12.500 days, wall time: 842.785 ms, max(u): (1.471e+00, 1.285e+00, 5.120e-03) m/s, next Δt: 20 minutes
[69.44%] i: 1000, t: 13.889 days, wall time: 919.861 ms, max(u): (1.419e+00, 1.176e+00, 4.376e-03) m/s, next Δt: 20 minutes
[76.39%] i: 1100, t: 15.278 days, wall time: 1.041 seconds, max(u): (1.597e+00, 1.256e+00, 3.696e-03) m/s, next Δt: 20 minutes
[83.33%] i: 1200, t: 16.667 days, wall time: 839.403 ms, max(u): (1.341e+00, 1.234e+00, 4.575e-03) m/s, next Δt: 20 minutes
[90.28%] i: 1300, t: 18.056 days, wall time: 831.399 ms, max(u): (1.291e+00, 1.378e+00, 2.820e-03) m/s, next Δt: 20 minutes
[97.22%] i: 1400, t: 19.444 days, wall time: 848.055 ms, max(u): (1.327e+00, 1.337e+00, 2.811e-03) m/s, next Δt: 20 minutes
[ Info: Simulation is stopping after running for 37.590 seconds.
[ Info: Simulation time 20 days equals or exceeds stop time 20 days.
[ Info: Simulation completed in 37.636 seconds
Visualization
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 $b$ and vertical vorticity $ζ$ at the surface, as well as the zonally-averaged zonal and meridional velocities $U$ and $V$ in the $(y, z)$ plane. First we load the 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.009373693726956844 -0.008098606020212173 -0.006873562000691891 -0.0056161051616072655 -0.004371653776615858 -0.0031541376374661922 -0.0018692617304623127 -0.0006117848097346723; -0.009386862628161907 -0.008134459145367146 -0.006854533217847347 -0.005626159254461527 -0.004373673349618912 -0.0031256743241101503 -0.0018897149711847305 -0.0006189537816680968; -0.009371398016810417 -0.008136555552482605 -0.006850988604128361 -0.005633965600281954 -0.004390474408864975 -0.0031390965450555086 -0.0018776304787024856 -0.0006329037714749575; -0.00937557127326727 -0.008122296072542667 -0.0068559893406927586 -0.005613374058157206 -0.004358673468232155 -0.0031190472654998302 -0.0018990880344063044 -0.000604456290602684; -0.009362660348415375 -0.008119618520140648 -0.006866539362818003 -0.005625402554869652 -0.0044029285199940205 -0.003105524927377701 -0.0018773219780996442 -0.0006168512045405805; -0.009378835558891296 -0.008125558495521545 -0.00684848427772522 -0.0056159961968660355 -0.00436010118573904 -0.00311826984398067 -0.0018437275430187583 -0.0006182187935337424; -0.009358277544379234 -0.008114620111882687 -0.006870206445455551 -0.0056203752756118774 -0.004371799994260073 -0.003131298581138253 -0.001888593309558928 -0.0006126259686425328; -0.009369920007884502 -0.00813291035592556 -0.006845311727374792 -0.005620083771646023 -0.004406928550451994 -0.003137404564768076 -0.0018858752446249127 -0.0006167318788357079; -0.009374761022627354 -0.008140740916132927 -0.006877920590341091 -0.005631486419588327 -0.004365687258541584 -0.0031087305396795273 -0.0018585214857012033 -0.0006290943128988147; -0.009373965673148632 -0.008140383288264275 -0.0068816994316875935 -0.005627751816064119 -0.0043582553043961525 -0.003134871833026409 -0.0018839569529518485 -0.0006621009670197964; -0.009435455314815044 -0.008134362287819386 -0.006861598696559668 -0.005617892369627953 -0.004376002587378025 -0.003115768311545253 -0.0018791657639667392 -0.0006194392335601151; -0.009404665790498257 -0.008116408251225948 -0.006854224484413862 -0.005618277937173843 -0.004404864739626646 -0.0031375959515571594 -0.0018738653743639588 -0.0006068337243050337; -0.009387591853737831 -0.00810569804161787 -0.006874336861073971 -0.005603729747235775 -0.004361129831522703 -0.0031116949394345284 -0.0018913225503638387 -0.0006397492252290249; -0.009380854666233063 -0.008121106773614883 -0.006890913937240839 -0.005610035732388496 -0.0043588788248598576 -0.0031099854968488216 -0.0018666110699996352 -0.0006441287696361542; -0.009365095756947994 -0.008147045038640499 -0.0068627651780843735 -0.005637490190565586 -0.00436483696103096 -0.00312750949524343 -0.0018543734913691878 -0.000614073418546468; -0.009366453625261784 -0.008137542754411697 -0.006891622673720121 -0.005604402162134647 -0.0043361796997487545 -0.0031129850540310144 -0.001865547732450068 -0.0006379625992849469; -0.009388178586959839 -0.0081177968531847 -0.006877831183373928 -0.00562796276062727 -0.0043618506751954556 -0.003116296138614416 -0.0018884784076362848 -0.0006182158249430358; -0.009379858151078224 -0.0081381406635046 -0.006883665919303894 -0.00564323365688324 -0.004387673921883106 -0.0031286957673728466 -0.0018726582638919353 -0.000641446269582957; -0.009363101795315742 -0.008126335218548775 -0.006869347766041756 -0.0056083910167217255 -0.0043652779422700405 -0.0030988790094852448 -0.0019031561678275466 -0.0006236854242160916; -0.009366431273519993 -0.008137797936797142 -0.006846864242106676 -0.0056292032822966576 -0.004386153072118759 -0.003122081747278571 -0.001866245293058455 -0.0006294718477874994; -0.009394514374434948 -0.008124042302370071 -0.00688911322504282 -0.005637291353195906 -0.004352377261966467 -0.003130368422716856 -0.0018727162387222052 -0.0006544884527102113; -0.009399520233273506 -0.0081134382635355 -0.006878588814288378 -0.005637849681079388 -0.0043650721199810505 -0.0031300268601626158 -0.0018589551327750087 -0.0006125365616753697; -0.007483987137675285 -0.00623729545623064 -0.004970815032720566 -0.0037393257953226566 -0.002505750395357609 -0.0012570226099342108 -9.346379556518514e-6 0.0012470358051359653; -0.00542523805052042 -0.0041806395165622234 -0.002929795766249299 -0.0016649016179144382 -0.0004086762492079288 0.0008141662110574543 0.0021141530014574528 0.0033211717382073402; -0.0033507337793707848 -0.002084780717268586 -0.0008548684418201447 0.00041193474316969514 0.001663030474446714 0.002903809305280447 0.00418754480779171 0.005395660176873207; -0.0012475264957174659 -1.1754656952689402e-5 0.0012463253224268556 0.002484834985807538 0.0037420678418129683 0.005002528429031372 0.006234708242118359 0.007522846572101116; 0.0006186122773215175 0.0018681746441870928 0.0031213609036058187 0.004365045577287674 0.005601604003459215 0.006873926613479853 0.008142880164086819 0.009395351633429527; 0.0006218919297680259 0.0018671528669074178 0.003145019058138132 0.004389484412968159 0.005612021312117577 0.006852028891444206 0.008119246922433376 0.009360314346849918; 0.0006430756184272468 0.0018866764148697257 0.0031020469032227993 0.0043753162026405334 0.005622196942567825 0.00688056880608201 0.00814918614923954 0.009371127933263779; 0.0006028130301274359 0.0018848218023777008 0.0031271805055439472 0.004353568889200687 0.005649688187986612 0.006875168066471815 0.008126472122967243 0.00937449187040329; 0.0006028561620041728 0.0018654277082532644 0.003094451269134879 0.004367747809737921 0.005632041487842798 0.006885922979563475 0.008145085535943508 0.009375941939651966; 0.0006085261702537537 0.0018853185465559363 0.0031141501385718584 0.004368348978459835 0.005656440742313862 0.006871630437672138 0.008116315118968487 0.00938541628420353; 0.0006205496028997004 0.0018865204183384776 0.00311265024356544 0.004384404979646206 0.005642211064696312 0.006895721890032291 0.00815179105848074 0.00935661792755127; 0.0006320195971056819 0.0018625101074576378 0.0031292159110307693 0.004345016088336706 0.00562835019081831 0.006873345002532005 0.008101584389805794 0.009401996619999409; 0.0006265567499212921 0.0018521146848797798 0.0031296364031732082 0.004364897031337023 0.00563524616882205 0.006867750082165003 0.008117998018860817 0.009406405501067638; 0.0006190070998854935 0.0018914623651653528 0.0031353835947811604 0.004364424850791693 0.005620671436190605 0.006893012672662735 0.008130027912557125 0.009370981715619564; 0.0006215891917236149 0.0018884632736444473 0.003138278378173709 0.004417513031512499 0.005594099406152964 0.006896191742271185 0.008113331161439419 0.009361766278743744; 0.0006527956575155258 0.0018975320272147655 0.003112254198640585 0.004362649284303188 0.005618505645543337 0.006861311849206686 0.008096783421933651 0.009388439357280731; 0.0006170784472487867 0.0018897855188697577 0.003137962194159627 0.004364623222500086 0.005630296654999256 0.006897061597555876 0.008108584210276604 0.009383193217217922; 0.0005998224369250238 0.0018726253183558583 0.003123372560366988 0.004372399300336838 0.005635806359350681 0.006883287336677313 0.00812196172773838 0.00936631765216589; 0.0006195336463861167 0.0018600632902234793 0.003105951240286231 0.004385001957416534 0.005610055755823851 0.006914208177477121 0.008144279941916466 0.009361697360873222; 0.0006032680394127965 0.0018881994765251875 0.00314209028147161 0.00438415864482522 0.005612269975244999 0.0068796975538134575 0.008137145079672337 0.009381057694554329; 0.0006533738342113793 0.0018700923537835479 0.0031128376722335815 0.004378338344395161 0.005619106814265251 0.006870719604194164 0.008119410835206509 0.009365593083202839; 0.0006291579338721931 0.0018481174483895302 0.00309008383192122 0.004377145320177078 0.005645364988595247 0.006869054399430752 0.008106213994324207 0.009372647851705551; 0.0006132190464995801 0.0018476187251508236 0.0031131806317716837 0.004350367467850447 0.005599306430667639 0.006867853458970785 0.008131074719130993 0.009380425326526165; 0.0006372373318299651 0.0018501480808481574 0.0031026836950331926 0.004374404437839985 0.005617949180305004 0.006893299985677004 0.008141743019223213 0.009372209198772907; 0.0006188727566041052 0.0018859871197491884 0.0031435233540832996 0.004379009362310171 0.0056202043779194355 0.006860880181193352 0.008104520849883556 0.00940529815852642; 0.000619645812548697 0.0018720864318311214 0.003121031913906336 0.004369039554148912 0.005634792614728212 0.006891835015267134 0.008146106265485287 0.009348972700536251])
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-29082/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-29082/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-29082/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.11
[9e8cae18] Oceananigans v0.104.2 `..`
[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.