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: 9.371 seconds, max(u): (0.000e+00, 0.000e+00, 0.000e+00) m/s, next Δt: 20 minutes
[ Info: ... simulation initialization complete (18.013 seconds)
[ Info: Executing initial time step...
[ Info: ... initial time step complete (4.510 seconds).
[06.94%] i: 100, t: 1.389 days, wall time: 21.551 seconds, max(u): (1.244e-01, 1.263e-01, 1.706e-03) m/s, next Δt: 20 minutes
[13.89%] i: 200, t: 2.778 days, wall time: 1.709 seconds, max(u): (2.199e-01, 1.999e-01, 2.039e-03) m/s, next Δt: 20 minutes
[20.83%] i: 300, t: 4.167 days, wall time: 2.800 seconds, max(u): (2.981e-01, 2.621e-01, 1.866e-03) m/s, next Δt: 20 minutes
[27.78%] i: 400, t: 5.556 days, wall time: 4.077 seconds, max(u): (3.836e-01, 4.001e-01, 1.956e-03) m/s, next Δt: 20 minutes
[34.72%] i: 500, t: 6.944 days, wall time: 1.264 seconds, max(u): (4.736e-01, 6.116e-01, 2.388e-03) m/s, next Δt: 20 minutes
[41.67%] i: 600, t: 8.333 days, wall time: 1.026 seconds, max(u): (7.126e-01, 9.156e-01, 3.225e-03) m/s, next Δt: 20 minutes
[48.61%] i: 700, t: 9.722 days, wall time: 1.025 seconds, max(u): (1.170e+00, 1.233e+00, 4.414e-03) m/s, next Δt: 20 minutes
[55.56%] i: 800, t: 11.111 days, wall time: 898.732 ms, max(u): (1.458e+00, 1.262e+00, 6.296e-03) m/s, next Δt: 20 minutes
[62.50%] i: 900, t: 12.500 days, wall time: 1.133 seconds, max(u): (1.264e+00, 1.175e+00, 4.661e-03) m/s, next Δt: 20 minutes
[69.44%] i: 1000, t: 13.889 days, wall time: 1.671 seconds, max(u): (1.232e+00, 1.343e+00, 6.121e-03) m/s, next Δt: 20 minutes
[76.39%] i: 1100, t: 15.278 days, wall time: 1.509 seconds, max(u): (1.197e+00, 1.260e+00, 3.125e-03) m/s, next Δt: 20 minutes
[83.33%] i: 1200, t: 16.667 days, wall time: 1.230 seconds, max(u): (1.468e+00, 1.266e+00, 3.106e-03) m/s, next Δt: 20 minutes
[90.28%] i: 1300, t: 18.056 days, wall time: 1.772 seconds, max(u): (1.513e+00, 1.373e+00, 3.094e-03) m/s, next Δt: 20 minutes
[97.22%] i: 1400, t: 19.444 days, wall time: 1.623 seconds, max(u): (1.537e+00, 1.457e+00, 2.837e-03) m/s, next Δt: 20 minutes
[ Info: Simulation is stopping after running for 46.317 seconds.
[ Info: Simulation time 20 days equals or exceeds stop time 20 days.
[ Info: Simulation completed in 46.424 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.009375235997140408 -0.00811192486435175 -0.00689184945076704 -0.005637611262500286 -0.004370285663753748 -0.0031166360713541508 -0.0019013759447261691 -0.0006241774535737932; -0.009377021342515945 -0.008126710541546345 -0.006897269748151302 -0.005583223421126604 -0.004373050760477781 -0.003112118924036622 -0.0018519897712394595 -0.0006350302137434483; -0.009384888224303722 -0.008128689602017403 -0.006841341499239206 -0.005627285689115524 -0.0043749394826591015 -0.0031025817152112722 -0.0018620709888637066 -0.0006331409094855189; -0.009375534951686859 -0.008111364208161831 -0.006878433283418417 -0.00565279321745038 -0.004350204486399889 -0.0031219515949487686 -0.0018732172902673483 -0.0006096221040934324; -0.009366188198328018 -0.00813918188214302 -0.006870608311146498 -0.005625865887850523 -0.004386974964290857 -0.0031068252865225077 -0.0018632888095453382 -0.000638713245280087; -0.00937227625399828 -0.008112931624054909 -0.006879803258925676 -0.0056379567831754684 -0.004386797081679106 -0.0031312063802033663 -0.0018676429754123092 -0.0006273716571740806; -0.009353190660476685 -0.008113563992083073 -0.00687221996486187 -0.0056014047004282475 -0.0043605235405266285 -0.0031289937905967236 -0.0018924153409898281 -0.0006224256940186024; -0.009391672909259796 -0.008121387101709843 -0.006846014875918627 -0.005617113318294287 -0.004395895637571812 -0.003120217937976122 -0.0018865966703742743 -0.0006149146356619895; -0.009358990006148815 -0.00812293030321598 -0.006884855218231678 -0.005623791366815567 -0.0043647014535963535 -0.003135467879474163 -0.0018773677293211222 -0.0006130700930953026; -0.009373990818858147 -0.008125358261168003 -0.006867343559861183 -0.0056161824613809586 -0.004369919188320637 -0.0031169112771749496 -0.0018638184992596507 -0.0006077504367567599; -0.009401408024132252 -0.008136020973324776 -0.006867650430649519 -0.005609733518213034 -0.004381394479423761 -0.0031221662648022175 -0.0019062922801822424 -0.0006326347938738763; -0.009373296052217484 -0.008143683895468712 -0.006875383201986551 -0.005592995323240757 -0.004374818876385689 -0.0031302417628467083 -0.0018730142619460821 -0.0006547613884322345; -0.009374154731631279 -0.008144269697368145 -0.006884158588945866 -0.005637274589389563 -0.0043848794884979725 -0.003134005004540086 -0.0018881210125982761 -0.000631016562692821; -0.009360553696751595 -0.008128282614052296 -0.006880699656903744 -0.005600098520517349 -0.004373607691377401 -0.0031430514063686132 -0.0019123409874737263 -0.0006167675019241869; -0.009369884617626667 -0.008126480504870415 -0.006891487166285515 -0.005611496977508068 -0.004381757695227861 -0.00312484847381711 -0.001864248188212514 -0.0006335237994790077; -0.009377380833029747 -0.008143550716340542 -0.006869255565106869 -0.00561570143327117 -0.004375815391540527 -0.003125069197267294 -0.0018622481729835272 -0.0006279026274569333; -0.009392758831381798 -0.008113470859825611 -0.006886160466820002 -0.00563042564317584 -0.004373256582766771 -0.0031142833177000284 -0.0018672056030482054 -0.0006175589514896274; -0.009361972101032734 -0.008130031637847424 -0.0068868608213961124 -0.005613674409687519 -0.004389950539916754 -0.003126595402136445 -0.0018794434145092964 -0.0006411268259398639; -0.009367265738546848 -0.008117719553411007 -0.006895359139889479 -0.0056320056319236755 -0.004386106505990028 -0.0031280608382076025 -0.0018780224490910769 -0.00062684336444363; -0.009359261021018028 -0.00812617875635624 -0.00685792975127697 -0.005625484045594931 -0.004354917909950018 -0.0031149734277278185 -0.0018810940673574805 -0.000629228597972542; -0.00935743935406208 -0.008125070482492447 -0.006887735333293676 -0.005656600464135408 -0.004389207344502211 -0.0031176009215414524 -0.00186338578350842 -0.0006226867553777993; -0.009386505000293255 -0.008118999190628529 -0.006873920559883118 -0.005643087439239025 -0.004371959716081619 -0.0031166153494268656 -0.001865257858298719 -0.0006132654380053282; -0.007518467027693987 -0.006256754510104656 -0.005012104753404856 -0.0037658619694411755 -0.002501204377040267 -0.0012308073928579688 2.7563105504668783e-6 0.001272544264793396; -0.005405817646533251 -0.004146842751652002 -0.002915803575888276 -0.001679415232501924 -0.00039878670941106975 0.0008460588287562132 0.002079383470118046 0.0033121060114353895; -0.003339701797813177 -0.0020785678643733263 -0.0008433181210421026 0.00040721523691900074 0.0016609436133876443 0.002906348090618849 0.004170320462435484 0.005383492913097143; -0.0012522997567430139 -4.448067556950264e-6 0.0012734616175293922 0.002492565428838134 0.00373278371989727 0.0050010476261377335 0.006254169624298811 0.0075123910792171955; 0.0006303091649897397 0.0018735735211521387 0.003109751734882593 0.00437902519479394 0.005624350626021624 0.006853856146335602 0.00812335405498743 0.009368356317281723; 0.0006209337734617293 0.0018964369082823396 0.0030933234374970198 0.004385674372315407 0.005635473877191544 0.006881017703562975 0.008117719553411007 0.009367161430418491; 0.0006285623530857265 0.0018749780720099807 0.0031119349878281355 0.004380548372864723 0.005613868124783039 0.006869707256555557 0.008136054500937462 0.009350734762847424; 0.0005968567566014826 0.0018690349534153938 0.0031022154726088047 0.004401927348226309 0.005610624793916941 0.006877702660858631 0.008138499222695827 0.009360717609524727; 0.0006225491524673998 0.00188490841537714 0.003135253209620714 0.004391293972730637 0.005623809061944485 0.006884084548801184 0.00813047681003809 0.00938584003597498; 0.0006280806264840066 0.0018825660226866603 0.0031156488694250584 0.004369994625449181 0.005633198656141758 0.0068550207652151585 0.008144090883433819 0.009382186457514763; 0.0006354746874421835 0.0018711359007284045 0.0031326718162745237 0.004401824902743101 0.005635466892272234 0.006863998249173164 0.00813811831176281 0.009384403936564922; 0.00061156100127846 0.0018824301660060883 0.003126903437077999 0.004377730656415224 0.005612218752503395 0.006904146634042263 0.008118100464344025 0.009347456507384777; 0.0006387713947333395 0.001876829075627029 0.003122229827567935 0.004382400307804346 0.005625233985483646 0.00688327057287097 0.00812057126313448 0.009381057694554329; 0.0006323278066702187 0.0018526638159528375 0.0031249746680259705 0.004385130945593119 0.005636267364025116 0.006861413363367319 0.00812328141182661 0.009392013773322105; 0.0006357469828799367 0.001890379935503006 0.0031583853997290134 0.00435087364166975 0.0056367372162640095 0.006875752937048674 0.008118907921016216 0.009367554448544979; 0.0006154570146463811 0.0018771664472296834 0.003130362369120121 0.004371595103293657 0.005619202274829149 0.00687048165127635 0.00810602493584156 0.009359845891594887; 0.000626249413471669 0.0019007653463631868 0.003130546072497964 0.004354455508291721 0.005629833787679672 0.006891061086207628 0.008142167702317238 0.00937369093298912; 0.0006191122811287642 0.0018779345555230975 0.0031199194490909576 0.004385385196655989 0.0056174322962760925 0.006881330162286758 0.008128187619149685 0.009356269612908363; 0.0006186654209159315 0.001899914932437241 0.003114563412964344 0.004373720847070217 0.0056296237744390965 0.006889224983751774 0.00811030063778162 0.009358111768960953; 0.0006307368166744709 0.0018962349276989698 0.0031018259469419718 0.0043886280618608 0.00561111094430089 0.006856165826320648 0.008104496635496616 0.00936251413077116; 0.0006211583968251944 0.0018484130268916488 0.0031156931072473526 0.004402745980769396 0.005630893167108297 0.006880922242999077 0.008132283575832844 0.009371495805680752; 0.0006408769986592233 0.0018996305298060179 0.003145522205159068 0.004408545792102814 0.005637818481773138 0.006885075941681862 0.008141531608998775 0.009386940859258175; 0.0006271499441936612 0.0018585611833259463 0.003102808492258191 0.004377047531306744 0.0056002177298069 0.006879081483930349 0.00812418945133686 0.009384021162986755; 0.0006095374119468033 0.0018834181828424335 0.0031078578904271126 0.0043729860335588455 0.005619900766760111 0.006878897547721863 0.008126514963805676 0.009396177716553211; 0.0006180992932058871 0.0018935676198452711 0.0031262030825018883 0.004367217421531677 0.005640086717903614 0.006862119305878878 0.00812746211886406 0.00937774870544672; 0.0006118855671957135 0.0018886416219174862 0.003148629330098629 0.004391435533761978 0.005639370065182447 0.006850803270936012 0.00812909472733736 0.009383853524923325])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-29665/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-29665/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-29665/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.