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: 2.0 MiBNow 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: 32.321 seconds, max(u): (0.000e+00, 0.000e+00, 0.000e+00) m/s, next Δt: 20 minutes
[ Info: ... simulation initialization complete (10.661 seconds)
[ Info: Executing initial time step...
[ Info: ... initial time step complete (3.964 seconds).
[06.94%] i: 100, t: 1.389 days, wall time: 12.152 seconds, max(u): (1.272e-01, 1.183e-01, 1.519e-03) m/s, next Δt: 20 minutes
[13.89%] i: 200, t: 2.778 days, wall time: 815.433 ms, max(u): (2.069e-01, 1.646e-01, 1.713e-03) m/s, next Δt: 20 minutes
[20.83%] i: 300, t: 4.167 days, wall time: 853.814 ms, max(u): (2.817e-01, 2.083e-01, 1.656e-03) m/s, next Δt: 20 minutes
[27.78%] i: 400, t: 5.556 days, wall time: 901.485 ms, max(u): (3.449e-01, 2.725e-01, 1.674e-03) m/s, next Δt: 20 minutes
[34.72%] i: 500, t: 6.944 days, wall time: 784.305 ms, max(u): (4.034e-01, 3.299e-01, 1.750e-03) m/s, next Δt: 20 minutes
[41.67%] i: 600, t: 8.333 days, wall time: 822.802 ms, max(u): (4.990e-01, 4.253e-01, 2.103e-03) m/s, next Δt: 20 minutes
[48.61%] i: 700, t: 9.722 days, wall time: 805.469 ms, max(u): (6.214e-01, 6.399e-01, 2.159e-03) m/s, next Δt: 20 minutes
[55.56%] i: 800, t: 11.111 days, wall time: 866.031 ms, max(u): (8.439e-01, 9.787e-01, 2.986e-03) m/s, next Δt: 20 minutes
[62.50%] i: 900, t: 12.500 days, wall time: 809.195 ms, max(u): (1.162e+00, 1.303e+00, 4.243e-03) m/s, next Δt: 20 minutes
[69.44%] i: 1000, t: 13.889 days, wall time: 843.622 ms, max(u): (1.451e+00, 1.310e+00, 4.309e-03) m/s, next Δt: 20 minutes
[76.39%] i: 1100, t: 15.278 days, wall time: 895.504 ms, max(u): (1.409e+00, 1.343e+00, 4.584e-03) m/s, next Δt: 20 minutes
[83.33%] i: 1200, t: 16.667 days, wall time: 819.766 ms, max(u): (1.311e+00, 1.472e+00, 4.504e-03) m/s, next Δt: 20 minutes
[90.28%] i: 1300, t: 18.056 days, wall time: 794.643 ms, max(u): (1.414e+00, 1.291e+00, 3.973e-03) m/s, next Δt: 20 minutes
[97.22%] i: 1400, t: 19.444 days, wall time: 845.288 ms, max(u): (1.416e+00, 1.179e+00, 3.755e-03) m/s, next Δt: 20 minutes
[ Info: Simulation is stopping after running for 29.143 seconds.
[ Info: Simulation time 20 days equals or exceeds stop time 20 days.
[ Info: Simulation completed in 29.180 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()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.009350836277008057 -0.008148481138050556 -0.0068885572254657745 -0.005620043724775314 -0.004361143335700035 -0.0031503031495958567 -0.0018633354920893908 -0.0006318832747638226; -0.009367200545966625 -0.008132304064929485 -0.006874625571072102 -0.0056384471245110035 -0.004370559472590685 -0.0031207362189888954 -0.0018649650737643242 -0.000625555869191885; -0.009392363019287586 -0.008132102899253368 -0.006873644422739744 -0.005627282429486513 -0.004367019049823284 -0.0031163778621703386 -0.0018565284553915262 -0.000630040536634624; -0.009360515512526035 -0.00813385285437107 -0.0068660094402730465 -0.005603002849966288 -0.004379598889499903 -0.0031161666847765446 -0.0018682156223803759 -0.0006118231103755534; -0.009397453628480434 -0.00816105306148529 -0.00687648169696331 -0.005616986192762852 -0.004385988228023052 -0.00312707107514143 -0.0018700208747759461 -0.0006083447369746864; -0.009382810443639755 -0.00810319185256958 -0.006875809282064438 -0.005648183170706034 -0.0043654851615428925 -0.0031154612079262733 -0.0018858907278627157 -0.0006230337312445045; -0.009388347156345844 -0.008137964643537998 -0.006872125435620546 -0.005625972058624029 -0.004349386785179377 -0.003120825160294771 -0.0018838883843272924 -0.0006291313911788166; -0.009353757835924625 -0.008108854293823242 -0.0068991235457360744 -0.005642585456371307 -0.004344810266047716 -0.003114827210083604 -0.0018362520495429635 -0.000613472075201571; -0.009361221455037594 -0.008103425614535809 -0.006855904124677181 -0.005643066484481096 -0.004390091635286808 -0.003111308440566063 -0.0018728150753304362 -0.0006383388536050916; -0.009378848597407341 -0.008142472244799137 -0.006867477670311928 -0.005612330976873636 -0.004340991377830505 -0.0031086946837604046 -0.0018652142025530338 -0.000621398095972836; -0.009395874105393887 -0.008118302561342716 -0.006847997661679983 -0.0056423889473080635 -0.004349029157310724 -0.0031214735936373472 -0.0018880346324294806 -0.0006231298320926726; -0.009381016716361046 -0.00812347698956728 -0.006883740425109863 -0.005620358977466822 -0.00438581733033061 -0.0031175899785012007 -0.0018951587844640017 -0.0005948653561063111; -0.009376473724842072 -0.008130112662911415 -0.006863177753984928 -0.005614027380943298 -0.00436791218817234 -0.0031151885632425547 -0.0018762934487313032 -0.0006104916683398187; -0.009369241073727608 -0.008118624798953533 -0.006863364018499851 -0.005645883735269308 -0.004381977021694183 -0.003086336422711611 -0.0018770082388073206 -0.0006274489569477737; -0.009388765320181847 -0.008135899901390076 -0.006868388503789902 -0.005642895586788654 -0.00437589967623353 -0.0031464998610317707 -0.001879603136330843 -0.0006171925342641771; -0.009367701597511768 -0.008114730939269066 -0.006890249438583851 -0.005624775309115648 -0.004371295683085918 -0.0031226808205246925 -0.0018725201953202486 -0.0006334466743282974; -0.009405364282429218 -0.00811870489269495 -0.006857398897409439 -0.0056635234504938126 -0.004358334466814995 -0.0031209418084472418 -0.001899904222227633 -0.000631842645816505; -0.009382542222738266 -0.008128199726343155 -0.006872630212455988 -0.00561985420063138 -0.004391056951135397 -0.0031160844955593348 -0.0018708441639319062 -0.0006366319139488041; -0.009373272769153118 -0.008136341348290443 -0.006835950538516045 -0.005605851300060749 -0.004413339775055647 -0.0031059379689395428 -0.0018674369202926755 -0.0006167740211822093; -0.009356578812003136 -0.008092761039733887 -0.006865420378744602 -0.005635654088109732 -0.00435701385140419 -0.003114215563982725 -0.0018674613675102592 -0.0006290691671893001; -0.009366222657263279 -0.008109874092042446 -0.00688145449385047 -0.005629145540297031 -0.004377664532512426 -0.003112517297267914 -0.0018921426963061094 -0.0006027946365065873; -0.009380919858813286 -0.008098188787698746 -0.006863737944513559 -0.005631303880363703 -0.004382694140076637 -0.0031271532643586397 -0.0018706700066104531 -0.0006222141091711819; -0.00749067310243845 -0.0062247696332633495 -0.00500487070530653 -0.003752992022782564 -0.0025098114274442196 -0.0012562540359795094 -5.103850526211318e-6 0.0012329204473644495; -0.005395224317908287 -0.0041780308820307255 -0.0029409765265882015 -0.0016404808266088367 -0.000397646042983979 0.0008556998218409717 0.0020612680818885565 0.0033531540539115667; -0.0033345522824674845 -0.0020664322655647993 -0.0008287822129204869 0.0004219897382427007 0.001670458703301847 0.002911404939368367 0.004180756863206625 0.005428378935903311; -0.0012791634071618319 9.255781151296105e-6 0.0012661590008065104 0.002491403603926301 0.00376697676256299 0.004980441648513079 0.006260776426643133 0.0074828690849244595; 0.0006095940480008721 0.001873429981060326 0.003108874661847949 0.004383117891848087 0.005610894411802292 0.006884552538394928 0.00811837986111641 0.009371146559715271; 0.0006242137169465423 0.0018460438586771488 0.00311819976195693 0.004373412579298019 0.005631246604025364 0.006875684019178152 0.008112265728414059 0.00939067080616951; 0.0006278369692154229 0.00188510084990412 0.0031216552015393972 0.004383408930152655 0.005608722567558289 0.006857092492282391 0.008105235174298286 0.009375149384140968; 0.0006202555960044265 0.0018816509982571006 0.0031192132737487555 0.004359052982181311 0.005621211603283882 0.006876248866319656 0.008175437338650227 0.009359169751405716; 0.0006365078152157366 0.0018584947101771832 0.0031444148626178503 0.004380891099572182 0.005611785221844912 0.006869015749543905 0.008130930364131927 0.009380843490362167; 0.0006101872422732413 0.0018881019204854965 0.003169191535562277 0.004350422881543636 0.005633524153381586 0.006868456490337849 0.008142028003931046 0.009351632557809353; 0.0006433036760427058 0.0018772901967167854 0.0031059742905199528 0.004361866973340511 0.005598554853349924 0.006883066613227129 0.008136214688420296 0.009360332041978836; 0.0006298546795733273 0.001860254444181919 0.0031174130272120237 0.004385391250252724 0.005649032536894083 0.006875913590192795 0.008118030615150928 0.009408644400537014; 0.0006270738085731864 0.0018766197608783841 0.0031340671703219414 0.004380003083497286 0.005607711151242256 0.006866972427815199 0.008146697655320168 0.00935191661119461; 0.0006295189377851784 0.0018629480618983507 0.0031244130805134773 0.004375657066702843 0.005659438669681549 0.006884139962494373 0.008135424926877022 0.009361420758068562; 0.000613624055404216 0.0018472423544153571 0.0031339111737906933 0.004379513207823038 0.005660144612193108 0.006866744253784418 0.008119416423141956 0.009389720857143402; 0.0006174884038046002 0.0018798013916239142 0.003139105625450611 0.004370252136141062 0.005646652076393366 0.006867650430649519 0.008138401433825493 0.009377519600093365; 0.0006029080832377076 0.0018845867598429322 0.0031134525779634714 0.004376575816422701 0.005624699871987104 0.0068747252225875854 0.008132348768413067 0.009348737075924873; 0.000621242739725858 0.0018867958569899201 0.0031349454075098038 0.004390019923448563 0.005634470377117395 0.006891743279993534 0.008150984533131123 0.009375999681651592; 0.000587008660659194 0.0018881328869611025 0.003123799804598093 0.004373593721538782 0.0056287129409611225 0.006862512324005365 0.008132797665894032 0.009370855987071991; 0.0006128349923528731 0.0018672043224796653 0.0031158658675849438 0.004366334527730942 0.005614263471215963 0.006878041662275791 0.008116800338029861 0.009373528882861137; 0.0006046831258572638 0.0018759341910481453 0.003111864672973752 0.004393294919282198 0.005653207655996084 0.006857128813862801 0.008136298507452011 0.009395243600010872; 0.0006284356932155788 0.0018880186835303903 0.0031422865577042103 0.004384628031402826 0.005612373352050781 0.006871083285659552 0.008104177191853523 0.009380675852298737; 0.0006235620239749551 0.0018763893749564886 0.0031121689826250076 0.004375392571091652 0.0056337155401706696 0.006856657098978758 0.008148475550115108 0.009354440495371819; 0.0005986347678117454 0.00188135274220258 0.0031497294548898935 0.004373400472104549 0.0056360638700425625 0.0068744332529604435 0.008113601244986057 0.00938639510422945; 0.0006196447066031396 0.0018852980574592948 0.0031243518460541964 0.004363468382507563 0.0056115505285561085 0.006889620795845985 0.008120749145746231 0.009379726834595203; 0.0006728464504703879 0.0018493300303816795 0.0031382928136736155 0.004384526051580906 0.005603456404060125 0.006877034902572632 0.0081336610019207 0.00937532540410757])
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
endThis page was generated using Literate.jl.