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.Units

Grid

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.0

Model

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(order=5)
│   └── b: WENO(order=5)
└── 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 . We impose the front on top of a vertical buoyancy gradient 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

# 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⁻²]")

fig

Simulation

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
├── Elapsed wall time: 0 seconds
├── Wall time per 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 entries
└── Diagnostics: OrderedDict with no entries

We 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: 32.5 KiB

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: 29.023 seconds, max(u): (0.000e+00, 0.000e+00, 0.000e+00) m/s, next Δt: 20 minutes
[ Info:     ... simulation initialization complete (27.096 seconds)
[ Info: Executing initial time step...
[ Info:     ... initial time step complete (19.788 seconds).
[06.94%] i: 100, t: 1.389 days, wall time: 38.859 seconds, max(u): (1.324e-01, 1.136e-01, 1.568e-03) m/s, next Δt: 20 minutes
[13.89%] i: 200, t: 2.778 days, wall time: 934.224 ms, max(u): (2.206e-01, 1.872e-01, 1.827e-03) m/s, next Δt: 20 minutes
[20.83%] i: 300, t: 4.167 days, wall time: 772.783 ms, max(u): (3.034e-01, 2.874e-01, 1.844e-03) m/s, next Δt: 20 minutes
[27.78%] i: 400, t: 5.556 days, wall time: 694.011 ms, max(u): (4.051e-01, 3.553e-01, 1.918e-03) m/s, next Δt: 20 minutes
[34.72%] i: 500, t: 6.944 days, wall time: 750.838 ms, max(u): (4.948e-01, 5.535e-01, 2.249e-03) m/s, next Δt: 20 minutes
[41.67%] i: 600, t: 8.333 days, wall time: 752.492 ms, max(u): (6.060e-01, 8.870e-01, 3.130e-03) m/s, next Δt: 20 minutes
[48.61%] i: 700, t: 9.722 days, wall time: 764.291 ms, max(u): (9.355e-01, 1.164e+00, 4.192e-03) m/s, next Δt: 20 minutes
[55.56%] i: 800, t: 11.111 days, wall time: 754.707 ms, max(u): (1.224e+00, 1.190e+00, 4.538e-03) m/s, next Δt: 20 minutes
[62.50%] i: 900, t: 12.500 days, wall time: 758.658 ms, max(u): (1.300e+00, 9.758e-01, 4.753e-03) m/s, next Δt: 20 minutes
[69.44%] i: 1000, t: 13.889 days, wall time: 833.219 ms, max(u): (1.245e+00, 9.535e-01, 3.434e-03) m/s, next Δt: 20 minutes
[76.39%] i: 1100, t: 15.278 days, wall time: 714.770 ms, max(u): (1.243e+00, 9.549e-01, 4.373e-03) m/s, next Δt: 20 minutes
[83.33%] i: 1200, t: 16.667 days, wall time: 651.677 ms, max(u): (1.317e+00, 1.195e+00, 3.646e-03) m/s, next Δt: 20 minutes
[90.28%] i: 1300, t: 18.056 days, wall time: 734.535 ms, max(u): (1.225e+00, 1.441e+00, 4.476e-03) m/s, next Δt: 20 minutes
[97.22%] i: 1400, t: 19.444 days, wall time: 669.787 ms, max(u): (1.519e+00, 1.303e+00, 3.180e-03) m/s, next Δt: 20 minutes
[ Info: Simulation is stopping after running for 0 seconds.
[ Info: Simulation time 20 days equals or exceeds stop time 20 days.
[ Info: Simulation completed in 1.023 minutes

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 CairoMakie

Three-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.grid
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.0

We 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)
41

Now 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
49-element Vector{Float64}:
 -500.0
 -479.1666666666667
 -458.3333333333333
 -437.5
 -416.6666666666667
 -395.8333333333333
 -375.0
 -354.1666666666667
 -333.3333333333333
 -312.5
 -291.6666666666667
 -270.8333333333333
 -250.0
 -229.16666666666666
 -208.33333333333334
 -187.5
 -166.66666666666666
 -145.83333333333334
 -125.0
 -104.16666666666667
  -83.33333333333333
  -62.5
  -41.666666666666664
  -20.833333333333332
    0.0
   20.833333333333332
   41.666666666666664
   62.5
   83.33333333333333
  104.16666666666667
  125.0
  145.83333333333334
  166.66666666666666
  187.5
  208.33333333333334
  229.16666666666666
  250.0
  270.8333333333333
  291.6666666666667
  312.5
  333.3333333333333
  354.1666666666667
  375.0
  395.8333333333333
  416.6666666666667
  437.5
  458.3333333333333
  479.1666666666667
  500.0

Next, 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!.

set_theme!(Theme(fontsize=24))

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.00936653558164835 -0.008112313225865364 -0.00686991261318326 -0.005612990353256464 -0.004373672418296337 -0.0031340920832008123 -0.0018624691292643547 -0.0006291430909186602; -0.009384457021951675 -0.008136201649904251 -0.006878307089209557 -0.005645302589982748 -0.0043863654136657715 -0.0031320233829319477 -0.0018485559849068522 -0.0006215755129233003; -0.009400860406458378 -0.008100182749330997 -0.006870339624583721 -0.005629376508295536 -0.004367370158433914 -0.0031512549612671137 -0.001872265012934804 -0.0006235940963961184; -0.009387461468577385 -0.008138267323374748 -0.006861629895865917 -0.00559549406170845 -0.004349303897470236 -0.0031358820851892233 -0.001881214790046215 -0.000624238746240735; -0.009371952153742313 -0.008088621310889721 -0.006875704973936081 -0.005626761820167303 -0.004380237776786089 -0.0031174567993730307 -0.0018708201823756099 -0.0006265028496272862; -0.009395972825586796 -0.008118095807731152 -0.006879355292767286 -0.0056231399066746235 -0.004349862691015005 -0.0031596655026078224 -0.0018826831365004182 -0.0006479555740952492; -0.009385774843394756 -0.008145410567522049 -0.0068677980452775955 -0.005630092229694128 -0.004389206413179636 -0.0031162211671471596 -0.0018677609041333199 -0.0006202551885508001; -0.009378268383443356 -0.0081388670951128 -0.006869935896247625 -0.005615886766463518 -0.004345064517110586 -0.003118587890639901 -0.0018847656901925802 -0.0005986696924082935; -0.009378398768603802 -0.008105969056487083 -0.006868243217468262 -0.00562689732760191 -0.0043912483379244804 -0.0031369635835289955 -0.0018726061098277569 -0.0006326762377284467; -0.009382651187479496 -0.008100167848169804 -0.006901987362653017 -0.005615381523966789 -0.004376746714115143 -0.0031197310891002417 -0.0018602295313030481 -0.000638815457932651; -0.00938786007463932 -0.008133401162922382 -0.0068664816208183765 -0.005619090516120195 -0.004367317538708448 -0.003113269340246916 -0.0018486259505152702 -0.000614236865658313; -0.009353636763989925 -0.008127862587571144 -0.006838849280029535 -0.005620231851935387 -0.004382755607366562 -0.003138534491881728 -0.0018909359350800514 -0.0006220031646080315; -0.009377200156450272 -0.008107300847768784 -0.006857138127088547 -0.005634082481265068 -0.004370416048914194 -0.0031270694453269243 -0.0018846163293346763 -0.000636632670648396; -0.009369520470499992 -0.008124944753944874 -0.006884061265736818 -0.005630732513964176 -0.004371758084744215 -0.0031171725131571293 -0.0018961651949211955 -0.0006092615658417344; -0.00936267152428627 -0.008122681640088558 -0.00690710823982954 -0.005634574219584465 -0.0043737213127315044 -0.003128479467704892 -0.0018788250163197517 -0.0006147946696728468; -0.009380357339978218 -0.008103903383016586 -0.006885112728923559 -0.005627789068967104 -0.004363317973911762 -0.003145291469991207 -0.0018584320787340403 -0.000619639118667692; -0.009390206076204777 -0.008117393590509892 -0.006853872444480658 -0.005613244604319334 -0.00437620934098959 -0.0031097782775759697 -0.0018720021471381187 -0.0006311176693998277; -0.009365624748170376 -0.008102347142994404 -0.006878937128931284 -0.005639639217406511 -0.004380390048027039 -0.0031284866854548454 -0.001868824940174818 -0.0006005216273479164; -0.00936424545943737 -0.00810347218066454 -0.006892349570989609 -0.00562502583488822 -0.004383375868201256 -0.0031264247372746468 -0.0018875249661505222 -0.0006485190242528915; -0.00939321331679821 -0.008118359372019768 -0.006883792579174042 -0.005616253707557917 -0.004382846411317587 -0.0031470651738345623 -0.0018489372450858355 -0.00061003677546978; -0.009360523894429207 -0.008134604431688786 -0.00687041413038969 -0.00562152499333024 -0.004352772142738104 -0.0030944363679736853 -0.0018810364417731762 -0.0006040629814378917; -0.009372216649353504 -0.008112729527056217 -0.006879628635942936 -0.005626433528959751 -0.0043667880818247795 -0.003121797926723957 -0.001868621213361621 -0.0006380933336913586; -0.007517815567553043 -0.0062231034971773624 -0.004989831242710352 -0.0037446131464093924 -0.0024849141482263803 -0.0012626455863937736 2.054678560625689e-7 0.0012648997362703085; -0.005413376726210117 -0.004160670097917318 -0.0029231926891952753 -0.0016646673902869225 -0.00040197346243076026 0.0008383152890019119 0.002056925091892481 0.0033228560350835323; -0.003328420454636216 -0.00207346398383379 -0.0008347588591277599 0.00041874434100463986 0.001664409413933754 0.0029251561500132084 0.004159871488809586 0.00539779057726264; -0.001235905452631414 2.270573895657435e-5 0.0012469172943383455 0.0024925312027335167 0.003746259957551956 0.005007854197174311 0.00623122975230217 0.007472809869796038; 0.0006302915280684829 0.0018749994924291968 0.0031026119831949472 0.004352365620434284 0.005654789041727781 0.006870476994663477 0.008149066008627415 0.009384733624756336; 0.0006453644018620253 0.0018798384116962552 0.0031260489486157894 0.004383433610200882 0.005627370439469814 0.006886512041091919 0.008116891607642174 0.009364870376884937; 0.0006334411446005106 0.0018486822955310345 0.0031477317679673433 0.004378946963697672 0.005641408264636993 0.00688225869089365 0.00814303569495678 0.00939034204930067; 0.0006258822395466268 0.0018938735593110323 0.003133761463686824 0.0043814824894070625 0.0056281983852386475 0.00687234103679657 0.008146960288286209 0.009389669634401798; 0.0006258225184865296 0.001873476430773735 0.0031178328208625317 0.004390884656459093 0.005632100626826286 0.006884823087602854 0.008135701529681683 0.009355691261589527; 0.0006068933289498091 0.0018676006002351642 0.0031262459233403206 0.0043734777718782425 0.005617150571197271 0.0068741911090910435 0.008124408312141895 0.009392224252223969; 0.0006231061415746808 0.001886315643787384 0.003116650739684701 0.004357154946774244 0.005614285357296467 0.006869203876703978 0.008122292347252369 0.009379883296787739; 0.0006241474184207618 0.0018683255184441805 0.0031337663531303406 0.004338475875556469 0.005638202652335167 0.006879940629005432 0.008122067898511887 0.009374825283885002; 0.0006074754055589437 0.0018508147913962603 0.0031081533525139093 0.004390641115605831 0.0056064315140247345 0.006863922346383333 0.008113227784633636 0.009374937042593956; 0.0006350050680339336 0.0018810504116117954 0.0031091361306607723 0.004384693689644337 0.0056120241060853004 0.006861833855509758 0.00813040230423212 0.009399698115885258; 0.0006028927164152265 0.0018722880631685257 0.003140214830636978 0.004375623073428869 0.005616519600152969 0.006879123859107494 0.008157224394381046 0.009349951520562172; 0.0005996755207888782 0.0018821748672053218 0.0031485718209296465 0.004381998907774687 0.005654904991388321 0.006880976725369692 0.008100584149360657 0.00936618447303772; 0.0006219056667760015 0.001841568504460156 0.0031315714586526155 0.004381267819553614 0.005623334553092718 0.006879906170070171 0.008127965033054352 0.009389129467308521; 0.0006037600105628371 0.0018523146864026785 0.0031394208781421185 0.004389877896755934 0.005624350160360336 0.006844925694167614 0.008115272037684917 0.009371084161102772; 0.0006120213074609637 0.0018937416607514024 0.0031377202831208706 0.004380404949188232 0.005619958508759737 0.006870955228805542 0.008142153732478619 0.009379632771015167; 0.0006102135521359742 0.0018871421925723553 0.0031242442782968283 0.004373392555862665 0.00561214704066515 0.0068842568434774876 0.008123001083731651 0.009385926648974419; 0.0006213795859366655 0.0018774475902318954 0.003128773532807827 0.004373460542410612 0.0056172870099544525 0.006856963504105806 0.008116449229419231 0.009360415861010551; 0.0006216756883077323 0.001888929051347077 0.003146365052089095 0.004401199985295534 0.005633274558931589 0.006886230781674385 0.008119790814816952 0.009387164376676083; 0.0006363036227412522 0.0018805282888934016 0.0031409801449626684 0.004407856147736311 0.005629687570035458 0.006852270103991032 0.008098848164081573 0.009376202709972858; 0.0006195693276822567 0.0018888147315010428 0.0031535779125988483 0.004372844938188791 0.005630128551274538 0.00687171146273613 0.00813211128115654 0.009360571391880512; 0.0006263646646402776 0.0018733496544882655 0.003116902429610491 0.004369563888758421 0.005650979932397604 0.006886712741106749 0.008136185817420483 0.009352924302220345; 0.0006313991034403443 0.001892293686978519 0.0031200870871543884 0.004387089051306248 0.005624009296298027 0.006841564085334539 0.0081389294937253 0.009371506050229073])

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
end


This page was generated using Literate.jl.