Single-column ocean simulation forced by JRA55 re-analysis
In this example, we simulate the evolution of an ocean water column forced by an atmosphere derived from the JRA55 re-analysis. The simulated column is located at ocean station Papa (144.9ᵒ W and 50.1ᵒ N)
Install dependencies
First let's make sure we have all required packages installed.
using Pkg
pkg"add Oceananigans, ClimaOcean, CairoMakie"
using ClimaOcean
using ClimaOcean.ECCO
using Oceananigans
using Oceananigans.Units
using Oceananigans.BuoyancyFormulations: buoyancy_frequency
using Oceananigans.Units: Time
using Dates
using Printf
WARNING: using Dates.Time in module ##225 conflicts with an existing identifier.
Construct the grid
First, we construct a single-column grid with 2 meter spacing located at ocean station Papa.
Ocean station papa location
location_name = "ocean_station_papa"
λ★, φ★ = 35.1, 50.1
grid = RectilinearGrid(size = 200,
x = λ★,
y = φ★,
z = (-400, 0),
topology = (Flat, Flat, Bounded))
1×1×200 RectilinearGrid{Float64, Oceananigans.Grids.Flat, Oceananigans.Grids.Flat, Oceananigans.Grids.Bounded} on Oceananigans.Architectures.CPU with 0×0×3 halo
├── Flat x = 35.1
├── Flat y = 50.1
└── Bounded z ∈ [-400.0, 0.0] regularly spaced with Δz=2.0
An "ocean simulation"
Next, we use ClimaOcean's ocean_simulation constructor to build a realistic ocean simulation on the single-column grid,
ocean = ocean_simulation(grid; Δt=10minutes, coriolis=FPlane(latitude = φ★))
Simulation of HydrostaticFreeSurfaceModel{CPU, RectilinearGrid}(time = 0 seconds, iteration = 0)
├── Next time step: 10 minutes
├── Elapsed wall time: 0 seconds
├── Wall time per iteration: NaN days
├── Stop time: Inf 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
which wraps around the ocean model
ocean.model
HydrostaticFreeSurfaceModel{CPU, RectilinearGrid}(time = 0 seconds, iteration = 0)
├── grid: 1×1×200 RectilinearGrid{Float64, Oceananigans.Grids.Flat, Oceananigans.Grids.Flat, Oceananigans.Grids.Bounded} on Oceananigans.Architectures.CPU with 0×0×3 halo
├── timestepper: QuasiAdamsBashforth2TimeStepper
├── tracers: (T, S, e)
├── closure: CATKEVerticalDiffusivity{VerticallyImplicitTimeDiscretization}
├── buoyancy: SeawaterBuoyancy with g=9.80665 and BoussinesqEquationOfState{Float64} with ĝ = NegativeZDirection()
├── advection scheme:
│ ├── momentum: Nothing
│ ├── T: Nothing
│ ├── S: Nothing
│ └── e: Nothing
└── coriolis: Oceananigans.Coriolis.FPlane{Float64}
We set initial conditions from ECCO:
set!(ocean.model, T=ECCOMetadatum(:temperature),
S=ECCOMetadatum(:salinity))
[ Info: Downloading ECCO data: temperature in /storage5/buildkite-agent/.julia-3483/scratchspaces/0376089a-ecfe-4b0e-a64f-9c555d74d754/ECCO...
[ Info: Downloading (size: 98.9 MiB)...
[ Info: ... downloaded 10.0 MiB (10% complete, 879.252 ms)
[ Info: ... downloaded 20.2 MiB (20% complete, 1.262 seconds)
[ Info: ... downloaded 29.8 MiB (30% complete, 1.545 seconds)
[ Info: ... downloaded 39.8 MiB (40% complete, 1.820 seconds)
[ Info: ... downloaded 49.5 MiB (50% complete, 2.076 seconds)
[ Info: ... downloaded 59.5 MiB (60% complete, 2.360 seconds)
[ Info: ... downloaded 69.4 MiB (70% complete, 2.613 seconds)
[ Info: ... downloaded 79.4 MiB (80% complete, 2.883 seconds)
[ Info: ... downloaded 89.3 MiB (90% complete, 3.155 seconds)
[ Info: ... downloaded 98.9 MiB (100% complete, 3.435 seconds)
[ Info: Inpainting ClimaOcean.DataWrangling.ECCO.ECCO4Monthly temperature data from 1992-01-01T00:00:00...
[ Info: ... (1.846 minutes)
[ Info: Downloading ECCO data: salinity in /storage5/buildkite-agent/.julia-3483/scratchspaces/0376089a-ecfe-4b0e-a64f-9c555d74d754/ECCO...
[ Info: Inpainting ClimaOcean.DataWrangling.ECCO.ECCO4Monthly salinity data from 1992-01-01T00:00:00...
[ Info: ... (1.852 minutes)
A prescribed atmosphere based on JRA55 re-analysis
We build a PrescribedAtmosphere at the same location as the single-colunm grid which is based on the JRA55 reanalysis.
atmosphere = JRA55PrescribedAtmosphere(longitude = λ★,
latitude = φ★,
end_date = DateTime(1990, 1, 31), # Last day of the simulation
backend = InMemory())
2×2×1×241 PrescribedAtmosphere{Float32} on Oceananigans.Grids.LatitudeLongitudeGrid:
├── times: 241-element StepRangeLen{Float64, Base.TwicePrecision{Float64}, Base.TwicePrecision{Float64}, Int64}
├── surface_layer_height: 10.0
└── boundary_layer_height: 512.0
This builds a representation of the atmosphere on the small grid
atmosphere.grid
2×2×1 LatitudeLongitudeGrid{Float32, Oceananigans.Grids.Bounded, Oceananigans.Grids.Bounded, Oceananigans.Grids.Flat} on Oceananigans.Architectures.CPU with 2×2×0 halo and with precomputed metrics
├── longitude: Bounded λ ∈ [34.5938, 35.7188] variably spaced with min(Δλ)=0.5625, max(Δλ)=0.5625
├── latitude: Bounded φ ∈ [49.4227, 50.5459] variably spaced with min(Δφ)=0.561619, max(Δφ)=0.561623
└── z: Flat z
Let's take a look at the atmospheric state
ua = interior(atmosphere.velocities.u, 1, 1, 1, :)
va = interior(atmosphere.velocities.v, 1, 1, 1, :)
Ta = interior(atmosphere.tracers.T, 1, 1, 1, :)
qa = interior(atmosphere.tracers.q, 1, 1, 1, :)
t_days = atmosphere.times[1:length(ua)] / days
using CairoMakie
set_theme!(Theme(linewidth=3, fontsize=24))
fig = Figure(size=(800, 1000))
axu = Axis(fig[2, 1], xlabel="Days since Jan 1 1990", ylabel="Atmosphere \n velocity (m s⁻¹)")
axT = Axis(fig[3, 1], xlabel="Days since Jan 1 1990", ylabel="Atmosphere \n temperature (ᵒK)")
axq = Axis(fig[4, 1], xlabel="Days since Jan 1 1990", ylabel="Atmosphere \n specific humidity")
Label(fig[1, 1], "Atmospheric state over ocean station Papa", tellwidth=false)
lines!(axu, t_days, ua, label="Zonal velocity")
lines!(axu, t_days, va, label="Meridional velocity")
ylims!(axu, -6, 6)
axislegend(axu, framevisible=false, nbanks=2, position=:lb)
lines!(axT, t_days, Ta)
lines!(axq, t_days, qa)
current_figure()
We continue constructing a simulation.
radiation = Radiation()
coupled_model = OceanSeaIceModel(ocean; atmosphere, radiation)
simulation = Simulation(coupled_model, Δt=ocean.Δt, stop_time=30days)
wall_clock = Ref(time_ns())
function progress(sim)
msg = "Ocean Station Papa"
msg *= string(", iter: ", iteration(sim), ", time: ", prettytime(sim))
elapsed = 1e-9 * (time_ns() - wall_clock[])
msg *= string(", wall time: ", prettytime(elapsed))
wall_clock[] = time_ns()
u, v, w = sim.model.ocean.model.velocities
msg *= @sprintf(", max|u|: (%.2e, %.2e)", maximum(abs, u), maximum(abs, v))
T = sim.model.ocean.model.tracers.T
S = sim.model.ocean.model.tracers.S
e = sim.model.ocean.model.tracers.e
τx = first(sim.model.interfaces.net_fluxes.ocean_surface.u)
τy = first(sim.model.interfaces.net_fluxes.ocean_surface.v)
Q = first(sim.model.interfaces.net_fluxes.ocean_surface.Q)
u★ = sqrt(sqrt(τx^2 + τy^2))
Nz = size(T, 3)
msg *= @sprintf(", u★: %.2f m s⁻¹", u★)
msg *= @sprintf(", Q: %.2f W m⁻²", Q)
msg *= @sprintf(", T₀: %.2f ᵒC", first(interior(T, 1, 1, Nz)))
msg *= @sprintf(", extrema(T): (%.2f, %.2f) ᵒC", minimum(T), maximum(T))
msg *= @sprintf(", S₀: %.2f g/kg", first(interior(S, 1, 1, Nz)))
msg *= @sprintf(", e₀: %.2e m² s⁻²", first(interior(e, 1, 1, Nz)))
@info msg
return nothing
end
simulation.callbacks[:progress] = Callback(progress, IterationInterval(100))
Callback of progress on IterationInterval(100)
Build flux outputs
τx = simulation.model.interfaces.net_fluxes.ocean_surface.u
τy = simulation.model.interfaces.net_fluxes.ocean_surface.v
JT = simulation.model.interfaces.net_fluxes.ocean_surface.T
Js = simulation.model.interfaces.net_fluxes.ocean_surface.S
E = simulation.model.interfaces.atmosphere_ocean_interface.fluxes.water_vapor
Qc = simulation.model.interfaces.atmosphere_ocean_interface.fluxes.sensible_heat
Qv = simulation.model.interfaces.atmosphere_ocean_interface.fluxes.latent_heat
ρₒ = simulation.model.interfaces.ocean_properties.reference_density
cₚ = simulation.model.interfaces.ocean_properties.heat_capacity
Q = ρₒ * cₚ * JT
ρτx = ρₒ * τx
ρτy = ρₒ * τy
N² = buoyancy_frequency(ocean.model)
κc = ocean.model.diffusivity_fields.κc
fluxes = (; ρτx, ρτy, E, Js, Qv, Qc)
auxiliary_fields = (; N², κc)
fields = merge(ocean.model.velocities, ocean.model.tracers, auxiliary_fields)
(u = 1×1×200 Field{Oceananigans.Grids.Face, Oceananigans.Grids.Center, Oceananigans.Grids.Center} on Oceananigans.Grids.RectilinearGrid on Oceananigans.Architectures.CPU
├── grid: 1×1×200 RectilinearGrid{Float64, Oceananigans.Grids.Flat, Oceananigans.Grids.Flat, Oceananigans.Grids.Bounded} on Oceananigans.Architectures.CPU with 0×0×3 halo
├── boundary conditions: FieldBoundaryConditions
│ └── west: Nothing, east: Nothing, south: Nothing, north: Nothing, bottom: Flux, top: Flux, immersed: ZeroFlux
└── data: 1×1×206 OffsetArray(::Array{Float64, 3}, 1:1, 1:1, -2:203) with eltype Float64 with indices 1:1×1:1×-2:203
└── max=0.0, min=0.0, mean=0.0, v = 1×1×200 Field{Oceananigans.Grids.Center, Oceananigans.Grids.Face, Oceananigans.Grids.Center} on Oceananigans.Grids.RectilinearGrid on Oceananigans.Architectures.CPU
├── grid: 1×1×200 RectilinearGrid{Float64, Oceananigans.Grids.Flat, Oceananigans.Grids.Flat, Oceananigans.Grids.Bounded} on Oceananigans.Architectures.CPU with 0×0×3 halo
├── boundary conditions: FieldBoundaryConditions
│ └── west: Nothing, east: Nothing, south: Nothing, north: Nothing, bottom: Flux, top: Flux, immersed: ZeroFlux
└── data: 1×1×206 OffsetArray(::Array{Float64, 3}, 1:1, 1:1, -2:203) with eltype Float64 with indices 1:1×1:1×-2:203
└── max=0.0, min=0.0, mean=0.0, w = ZeroField{Int64}, T = 1×1×200 Field{Oceananigans.Grids.Center, Oceananigans.Grids.Center, Oceananigans.Grids.Center} on Oceananigans.Grids.RectilinearGrid on Oceananigans.Architectures.CPU
├── grid: 1×1×200 RectilinearGrid{Float64, Oceananigans.Grids.Flat, Oceananigans.Grids.Flat, Oceananigans.Grids.Bounded} on Oceananigans.Architectures.CPU with 0×0×3 halo
├── boundary conditions: FieldBoundaryConditions
│ └── west: Nothing, east: Nothing, south: Nothing, north: Nothing, bottom: ZeroFlux, top: Flux, immersed: ZeroFlux
└── data: 1×1×206 OffsetArray(::Array{Float64, 3}, 1:1, 1:1, -2:203) with eltype Float64 with indices 1:1×1:1×-2:203
└── max=17.5153, min=14.4472, mean=16.1953, S = 1×1×200 Field{Oceananigans.Grids.Center, Oceananigans.Grids.Center, Oceananigans.Grids.Center} on Oceananigans.Grids.RectilinearGrid on Oceananigans.Architectures.CPU
├── grid: 1×1×200 RectilinearGrid{Float64, Oceananigans.Grids.Flat, Oceananigans.Grids.Flat, Oceananigans.Grids.Bounded} on Oceananigans.Architectures.CPU with 0×0×3 halo
├── boundary conditions: FieldBoundaryConditions
│ └── west: Nothing, east: Nothing, south: Nothing, north: Nothing, bottom: ZeroFlux, top: Flux, immersed: ZeroFlux
└── data: 1×1×206 OffsetArray(::Array{Float64, 3}, 1:1, 1:1, -2:203) with eltype Float64 with indices 1:1×1:1×-2:203
└── max=39.1557, min=38.6286, mean=38.9938, e = 1×1×200 Field{Oceananigans.Grids.Center, Oceananigans.Grids.Center, Oceananigans.Grids.Center} on Oceananigans.Grids.RectilinearGrid on Oceananigans.Architectures.CPU
├── grid: 1×1×200 RectilinearGrid{Float64, Oceananigans.Grids.Flat, Oceananigans.Grids.Flat, Oceananigans.Grids.Bounded} on Oceananigans.Architectures.CPU with 0×0×3 halo
├── boundary conditions: FieldBoundaryConditions
│ └── west: Nothing, east: Nothing, south: Nothing, north: Nothing, bottom: ZeroFlux, top: Flux, immersed: ZeroFlux
└── data: 1×1×206 OffsetArray(::Array{Float64, 3}, 1:1, 1:1, -2:203) with eltype Float64 with indices 1:1×1:1×-2:203
└── max=0.0, min=0.0, mean=0.0, N² = KernelFunctionOperation at (Center, Center, Face)
├── grid: 1×1×200 RectilinearGrid{Float64, Oceananigans.Grids.Flat, Oceananigans.Grids.Flat, Oceananigans.Grids.Bounded} on Oceananigans.Architectures.CPU with 0×0×3 halo
├── kernel_function: ∂z_b (generic function with 4 methods)
└── arguments: ("SeawaterBuoyancy with g=9.80665 and BoussinesqEquationOfState{Float64}", "(T=1×1×200 Field{Oceananigans.Grids.Center, Oceananigans.Grids.Center, Oceananigans.Grids.Center} on Oceananigans.Grids.RectilinearGrid on Oceananigans.Architectures.CPU, S=1×1×200 Field{Oceananigans.Grids.Center, Oceananigans.Grids.Center, Oceananigans.Grids.Center} on Oceananigans.Grids.RectilinearGrid on Oceananigans.Architectures.CPU, e=1×1×200 Field{Oceananigans.Grids.Center, Oceananigans.Grids.Center, Oceananigans.Grids.Center} on Oceananigans.Grids.RectilinearGrid on Oceananigans.Architectures.CPU)"), κc = 1×1×201 Field{Oceananigans.Grids.Center, Oceananigans.Grids.Center, Oceananigans.Grids.Face} on Oceananigans.Grids.RectilinearGrid on Oceananigans.Architectures.CPU
├── grid: 1×1×200 RectilinearGrid{Float64, Oceananigans.Grids.Flat, Oceananigans.Grids.Flat, Oceananigans.Grids.Bounded} on Oceananigans.Architectures.CPU with 0×0×3 halo
├── boundary conditions: FieldBoundaryConditions
│ └── west: Nothing, east: Nothing, south: Nothing, north: Nothing, bottom: Nothing, top: Nothing, immersed: ZeroFlux
└── data: 1×1×207 OffsetArray(::Array{Float64, 3}, 1:1, 1:1, -2:204) with eltype Float64 with indices 1:1×1:1×-2:204
└── max=7.16294e-5, min=0.0, mean=7.70266e-7)
Slice fields at the surface
outputs = merge(fields, fluxes)
filename = "single_column_omip_$(location_name)"
simulation.output_writers[:jld2] = JLD2Writer(ocean.model, outputs; filename,
schedule = TimeInterval(3hours),
overwrite_existing = true)
run!(simulation)
filename *= ".jld2"
u = FieldTimeSeries(filename, "u")
v = FieldTimeSeries(filename, "v")
T = FieldTimeSeries(filename, "T")
S = FieldTimeSeries(filename, "S")
e = FieldTimeSeries(filename, "e")
N² = FieldTimeSeries(filename, "N²")
κ = FieldTimeSeries(filename, "κc")
Qv = FieldTimeSeries(filename, "Qv")
Qc = FieldTimeSeries(filename, "Qc")
Js = FieldTimeSeries(filename, "Js")
Ev = FieldTimeSeries(filename, "E")
ρτx = FieldTimeSeries(filename, "ρτx")
ρτy = FieldTimeSeries(filename, "ρτy")
Nz = size(T, 3)
times = Qc.times
ua = atmosphere.velocities.u
va = atmosphere.velocities.v
Ta = atmosphere.tracers.T
qa = atmosphere.tracers.q
Qlw = atmosphere.downwelling_radiation.longwave
Qsw = atmosphere.downwelling_radiation.shortwave
Pr = atmosphere.freshwater_flux.rain
Ps = atmosphere.freshwater_flux.snow
Nt = length(times)
uat = zeros(Nt)
vat = zeros(Nt)
Tat = zeros(Nt)
qat = zeros(Nt)
Qswt = zeros(Nt)
Qlwt = zeros(Nt)
Pt = zeros(Nt)
for n = 1:Nt
t = times[n]
uat[n] = ua[1, 1, 1, Time(t)]
vat[n] = va[1, 1, 1, Time(t)]
Tat[n] = Ta[1, 1, 1, Time(t)]
qat[n] = qa[1, 1, 1, Time(t)]
Qswt[n] = Qsw[1, 1, 1, Time(t)]
Qlwt[n] = Qlw[1, 1, 1, Time(t)]
Pt[n] = Pr[1, 1, 1, Time(t)] + Ps[1, 1, 1, Time(t)]
end
fig = Figure(size=(1800, 1800))
axτ = Axis(fig[1, 1:3], xlabel="Days since Oct 1 1992", ylabel="Wind stress (N m⁻²)")
axQ = Axis(fig[1, 4:6], xlabel="Days since Oct 1 1992", ylabel="Heat flux (W m⁻²)")
axu = Axis(fig[2, 1:3], xlabel="Days since Oct 1 1992", ylabel="Velocities (m s⁻¹)")
axT = Axis(fig[2, 4:6], xlabel="Days since Oct 1 1992", ylabel="Surface temperature (ᵒC)")
axF = Axis(fig[3, 1:3], xlabel="Days since Oct 1 1992", ylabel="Freshwater volume flux (m s⁻¹)")
axS = Axis(fig[3, 4:6], xlabel="Days since Oct 1 1992", ylabel="Surface salinity (g kg⁻¹)")
axuz = Axis(fig[4:5, 1:2], xlabel="Velocities (m s⁻¹)", ylabel="z (m)")
axTz = Axis(fig[4:5, 3:4], xlabel="Temperature (ᵒC)", ylabel="z (m)")
axSz = Axis(fig[4:5, 5:6], xlabel="Salinity (g kg⁻¹)", ylabel="z (m)")
axNz = Axis(fig[6:7, 1:2], xlabel="Buoyancy frequency (s⁻²)", ylabel="z (m)")
axκz = Axis(fig[6:7, 3:4], xlabel="Eddy diffusivity (m² s⁻¹)", ylabel="z (m)", xscale=log10)
axez = Axis(fig[6:7, 5:6], xlabel="Turbulent kinetic energy (m² s⁻²)", ylabel="z (m)", xscale=log10)
title = @sprintf("Single-column simulation at %.2f, %.2f", φ★, λ★)
Label(fig[0, 1:6], title)
n = Observable(1)
times = (times .- times[1]) ./days
Nt = length(times)
tn = @lift times[$n]
colors = Makie.wong_colors()
ρₒ = coupled_model.interfaces.ocean_properties.reference_density
τx = interior(ρτx, 1, 1, 1, :) ./ ρₒ
τy = interior(ρτy, 1, 1, 1, :) ./ ρₒ
u★ = @. (τx^2 + τy^2)^(1/4)
lines!(axu, times, interior(u, 1, 1, Nz, :), color=colors[1], label="Zonal")
lines!(axu, times, interior(v, 1, 1, Nz, :), color=colors[2], label="Meridional")
lines!(axu, times, u★, color=colors[3], label="Ocean-side u★")
vlines!(axu, tn, linewidth=4, color=(:black, 0.5))
axislegend(axu)
lines!(axτ, times, interior(ρτx, 1, 1, 1, :), label="Zonal")
lines!(axτ, times, interior(ρτy, 1, 1, 1, :), label="Meridional")
vlines!(axτ, tn, linewidth=4, color=(:black, 0.5))
axislegend(axτ)
lines!(axT, times, Tat[1:Nt] .- 273.15, color=colors[1], linewidth=2, linestyle=:dash, label="Atmosphere temperature")
lines!(axT, times, interior(T, 1, 1, Nz, :), color=colors[2], linewidth=4, label="Ocean surface temperature")
vlines!(axT, tn, linewidth=4, color=(:black, 0.5))
axislegend(axT)
lines!(axQ, times, interior(Qv, 1, 1, 1, 1:Nt), color=colors[2], label="Sensible", linewidth=2)
lines!(axQ, times, interior(Qc, 1, 1, 1, 1:Nt), color=colors[3], label="Latent", linewidth=2)
lines!(axQ, times, - interior(Qsw, 1, 1, 1, 1:Nt), color=colors[4], label="Shortwave", linewidth=2)
lines!(axQ, times, - interior(Qlw, 1, 1, 1, 1:Nt), color=colors[5], label="Longwave", linewidth=2)
vlines!(axQ, tn, linewidth=4, color=(:black, 0.5))
axislegend(axQ)
lines!(axF, times, Pt[1:Nt], label="Prescribed freshwater flux")
lines!(axF, times, - interior(Ev, 1, 1, 1, 1:Nt), label="Evaporation")
vlines!(axF, tn, linewidth=4, color=(:black, 0.5))
axislegend(axF)
lines!(axS, times, interior(S, 1, 1, Nz, :))
vlines!(axS, tn, linewidth=4, color=(:black, 0.5))
zc = znodes(T)
zf = znodes(κ)
un = @lift interior(u[$n], 1, 1, :)
vn = @lift interior(v[$n], 1, 1, :)
Tn = @lift interior(T[$n], 1, 1, :)
Sn = @lift interior(S[$n], 1, 1, :)
κn = @lift interior(κ[$n], 1, 1, :)
en = @lift interior(e[$n], 1, 1, :)
N²n = @lift interior(N²[$n], 1, 1, :)
scatterlines!(axuz, un, zc, label="u")
scatterlines!(axuz, vn, zc, label="v")
scatterlines!(axTz, Tn, zc)
scatterlines!(axSz, Sn, zc)
scatterlines!(axez, en, zc)
scatterlines!(axNz, N²n, zf)
scatterlines!(axκz, κn, zf)
axislegend(axuz)
ulim = max(maximum(abs, u), maximum(abs, v))
xlims!(axuz, -ulim, ulim)
Tmax = maximum(interior(T))
Tmin = minimum(interior(T))
xlims!(axTz, Tmin - 0.1, Tmax + 0.1)
Nmax = maximum(interior(N²))
xlims!(axNz, -Nmax/10, Nmax * 1.05)
κmax = maximum(interior(κ))
xlims!(axκz, 1e-9, κmax * 1.1)
emax = maximum(interior(e))
xlims!(axez, 1e-11, emax * 1.1)
Smax = maximum(interior(S))
Smin = minimum(interior(S))
xlims!(axSz, Smin - 0.2, Smax + 0.2)
record(fig, "single_column_profiles.mp4", 1:Nt, framerate=24) do nn
@info "Drawing frame $nn of $Nt..."
n[] = nn
end
[ Info: Initializing simulation...
[ Info: Ocean Station Papa, iter: 0, time: 0 seconds, wall time: 14.174 seconds, max|u|: (0.00e+00, 0.00e+00), u★: 0.00 m s⁻¹, Q: 109.06 W m⁻², T₀: 16.97 ᵒC, extrema(T): (14.45, 17.52) ᵒC, S₀: 38.63 g/kg, e₀: 0.00e+00 m² s⁻²
[ Info: ... simulation initialization complete (6.849 seconds)
[ Info: Executing initial time step...
[ Info: ... initial time step complete (5.163 seconds).
[ Info: Ocean Station Papa, iter: 100, time: 16.667 hours, wall time: 9.377 seconds, max|u|: (4.46e-06, 2.85e-04), u★: 0.00 m s⁻¹, Q: 152.78 W m⁻², T₀: 16.96 ᵒC, extrema(T): (14.45, 17.52) ᵒC, S₀: 38.63 g/kg, e₀: 1.92e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 200, time: 1.389 days, wall time: 280.079 ms, max|u|: (1.72e-04, 5.71e-04), u★: 0.00 m s⁻¹, Q: 114.01 W m⁻², T₀: 16.94 ᵒC, extrema(T): (14.45, 17.52) ᵒC, S₀: 38.63 g/kg, e₀: 3.90e-05 m² s⁻²
[ Info: Ocean Station Papa, iter: 300, time: 2.083 days, wall time: 273.244 ms, max|u|: (1.31e-04, 2.46e-04), u★: 0.00 m s⁻¹, Q: 110.84 W m⁻², T₀: 16.92 ᵒC, extrema(T): (14.45, 17.52) ᵒC, S₀: 38.63 g/kg, e₀: 1.62e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 400, time: 2.778 days, wall time: 313.973 ms, max|u|: (4.64e-04, 3.56e-04), u★: 0.00 m s⁻¹, Q: 133.33 W m⁻², T₀: 16.91 ᵒC, extrema(T): (14.45, 17.52) ᵒC, S₀: 38.64 g/kg, e₀: 1.06e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 500, time: 3.472 days, wall time: 271.491 ms, max|u|: (5.86e-04, 7.28e-04), u★: 0.00 m s⁻¹, Q: 190.06 W m⁻², T₀: 16.90 ᵒC, extrema(T): (14.45, 17.52) ᵒC, S₀: 38.64 g/kg, e₀: 1.24e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 600, time: 4.167 days, wall time: 278.223 ms, max|u|: (2.48e-04, 2.87e-04), u★: 0.00 m s⁻¹, Q: 124.63 W m⁻², T₀: 16.88 ᵒC, extrema(T): (14.45, 17.52) ᵒC, S₀: 38.65 g/kg, e₀: 1.57e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 700, time: 4.861 days, wall time: 259.368 ms, max|u|: (1.44e-04, 4.03e-04), u★: 0.00 m s⁻¹, Q: 91.65 W m⁻², T₀: 16.87 ᵒC, extrema(T): (14.45, 17.52) ᵒC, S₀: 38.65 g/kg, e₀: 8.17e-05 m² s⁻²
[ Info: Ocean Station Papa, iter: 800, time: 5.556 days, wall time: 276.221 ms, max|u|: (3.01e-04, 3.15e-04), u★: 0.00 m s⁻¹, Q: 104.05 W m⁻², T₀: 16.86 ᵒC, extrema(T): (14.45, 17.52) ᵒC, S₀: 38.65 g/kg, e₀: 1.38e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 900, time: 6.250 days, wall time: 257.479 ms, max|u|: (5.69e-04, 2.72e-04), u★: 0.00 m s⁻¹, Q: 125.75 W m⁻², T₀: 16.85 ᵒC, extrema(T): (14.45, 17.52) ᵒC, S₀: 38.65 g/kg, e₀: 1.25e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 1000, time: 6.944 days, wall time: 281.825 ms, max|u|: (4.31e-04, 7.89e-06), u★: 0.00 m s⁻¹, Q: 94.63 W m⁻², T₀: 16.84 ᵒC, extrema(T): (14.45, 17.52) ᵒC, S₀: 38.65 g/kg, e₀: 1.39e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 1100, time: 7.639 days, wall time: 270.603 ms, max|u|: (7.30e-04, 2.59e-04), u★: 0.00 m s⁻¹, Q: 99.04 W m⁻², T₀: 16.83 ᵒC, extrema(T): (14.45, 17.52) ᵒC, S₀: 38.66 g/kg, e₀: 1.40e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 1200, time: 8.333 days, wall time: 268.412 ms, max|u|: (5.17e-04, 1.55e-04), u★: 0.00 m s⁻¹, Q: 168.81 W m⁻², T₀: 16.82 ᵒC, extrema(T): (14.45, 17.52) ᵒC, S₀: 38.66 g/kg, e₀: 9.51e-05 m² s⁻²
[ Info: Ocean Station Papa, iter: 1300, time: 9.028 days, wall time: 272.659 ms, max|u|: (4.57e-04, 3.44e-04), u★: 0.00 m s⁻¹, Q: 103.74 W m⁻², T₀: 16.81 ᵒC, extrema(T): (14.45, 17.52) ᵒC, S₀: 38.66 g/kg, e₀: 2.03e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 1400, time: 9.722 days, wall time: 267.356 ms, max|u|: (2.87e-04, 6.50e-04), u★: 0.00 m s⁻¹, Q: 94.98 W m⁻², T₀: 16.80 ᵒC, extrema(T): (14.45, 17.52) ᵒC, S₀: 38.67 g/kg, e₀: 1.53e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 1500, time: 10.417 days, wall time: 267.957 ms, max|u|: (1.10e-03, 1.31e-03), u★: 0.00 m s⁻¹, Q: 128.53 W m⁻², T₀: 16.79 ᵒC, extrema(T): (14.45, 17.52) ᵒC, S₀: 38.67 g/kg, e₀: 1.20e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 1600, time: 11.111 days, wall time: 252.539 ms, max|u|: (5.02e-04, 6.90e-04), u★: 0.00 m s⁻¹, Q: 124.38 W m⁻², T₀: 16.78 ᵒC, extrema(T): (14.45, 17.52) ᵒC, S₀: 38.68 g/kg, e₀: 1.98e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 1700, time: 11.806 days, wall time: 289.121 ms, max|u|: (1.30e-04, 7.40e-04), u★: 0.00 m s⁻¹, Q: 75.35 W m⁻², T₀: 16.77 ᵒC, extrema(T): (14.45, 17.52) ᵒC, S₀: 38.68 g/kg, e₀: 1.03e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 1800, time: 12.500 days, wall time: 271.825 ms, max|u|: (4.57e-04, 9.19e-04), u★: 0.00 m s⁻¹, Q: 162.28 W m⁻², T₀: 16.76 ᵒC, extrema(T): (14.45, 17.52) ᵒC, S₀: 38.68 g/kg, e₀: 1.22e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 1900, time: 13.194 days, wall time: 280.883 ms, max|u|: (4.46e-04, 9.16e-04), u★: 0.00 m s⁻¹, Q: 176.26 W m⁻², T₀: 16.74 ᵒC, extrema(T): (14.45, 17.52) ᵒC, S₀: 38.68 g/kg, e₀: 2.19e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 2000, time: 13.889 days, wall time: 278.436 ms, max|u|: (1.08e-04, 7.04e-04), u★: 0.00 m s⁻¹, Q: 199.46 W m⁻², T₀: 16.73 ᵒC, extrema(T): (14.45, 17.51) ᵒC, S₀: 38.69 g/kg, e₀: 2.74e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 2100, time: 14.583 days, wall time: 271.697 ms, max|u|: (1.01e-04, 7.04e-04), u★: 0.00 m s⁻¹, Q: 162.37 W m⁻², T₀: 16.72 ᵒC, extrema(T): (14.45, 17.49) ᵒC, S₀: 38.69 g/kg, e₀: 2.28e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 2200, time: 15.278 days, wall time: 296.706 ms, max|u|: (1.62e-04, 5.49e-04), u★: 0.00 m s⁻¹, Q: 213.00 W m⁻², T₀: 16.70 ᵒC, extrema(T): (14.45, 17.47) ᵒC, S₀: 38.70 g/kg, e₀: 1.84e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 2300, time: 15.972 days, wall time: 281.624 ms, max|u|: (3.69e-04, 4.12e-04), u★: 0.00 m s⁻¹, Q: 151.94 W m⁻², T₀: 16.68 ᵒC, extrema(T): (14.45, 17.47) ᵒC, S₀: 38.70 g/kg, e₀: 2.43e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 2400, time: 16.667 days, wall time: 316.754 ms, max|u|: (4.43e-04, 2.44e-04), u★: 0.00 m s⁻¹, Q: 124.06 W m⁻², T₀: 16.68 ᵒC, extrema(T): (14.45, 17.45) ᵒC, S₀: 38.70 g/kg, e₀: 2.00e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 2500, time: 17.361 days, wall time: 274.638 ms, max|u|: (4.73e-04, 1.11e-05), u★: 0.00 m s⁻¹, Q: 171.43 W m⁻², T₀: 16.67 ᵒC, extrema(T): (14.45, 17.45) ᵒC, S₀: 38.71 g/kg, e₀: 7.94e-05 m² s⁻²
[ Info: Ocean Station Papa, iter: 2600, time: 18.056 days, wall time: 290.541 ms, max|u|: (4.06e-04, 1.82e-04), u★: 0.00 m s⁻¹, Q: 161.80 W m⁻², T₀: 16.66 ᵒC, extrema(T): (14.45, 17.43) ᵒC, S₀: 38.71 g/kg, e₀: 2.45e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 2700, time: 18.750 days, wall time: 307.873 ms, max|u|: (2.73e-04, 3.10e-04), u★: 0.00 m s⁻¹, Q: 183.52 W m⁻², T₀: 16.65 ᵒC, extrema(T): (14.45, 17.43) ᵒC, S₀: 38.71 g/kg, e₀: 2.57e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 2800, time: 19.444 days, wall time: 290.057 ms, max|u|: (1.14e-04, 3.88e-04), u★: 0.00 m s⁻¹, Q: 165.79 W m⁻², T₀: 16.63 ᵒC, extrema(T): (14.45, 17.40) ᵒC, S₀: 38.72 g/kg, e₀: 6.04e-05 m² s⁻²
[ Info: Ocean Station Papa, iter: 2900, time: 20.139 days, wall time: 294.798 ms, max|u|: (7.30e-05, 3.62e-04), u★: 0.00 m s⁻¹, Q: 136.51 W m⁻², T₀: 16.62 ᵒC, extrema(T): (14.45, 17.40) ᵒC, S₀: 38.72 g/kg, e₀: 2.26e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 3000, time: 20.833 days, wall time: 277.588 ms, max|u|: (2.22e-04, 2.70e-04), u★: 0.00 m s⁻¹, Q: 164.66 W m⁻², T₀: 16.62 ᵒC, extrema(T): (14.45, 17.40) ᵒC, S₀: 38.72 g/kg, e₀: 2.59e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 3100, time: 21.528 days, wall time: 284.798 ms, max|u|: (3.14e-04, 1.50e-04), u★: 0.00 m s⁻¹, Q: 177.79 W m⁻², T₀: 16.61 ᵒC, extrema(T): (14.45, 17.38) ᵒC, S₀: 38.72 g/kg, e₀: 1.57e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 3200, time: 22.222 days, wall time: 289.990 ms, max|u|: (3.19e-04, 4.43e-06), u★: 0.00 m s⁻¹, Q: 133.06 W m⁻², T₀: 16.60 ᵒC, extrema(T): (14.45, 17.38) ᵒC, S₀: 38.73 g/kg, e₀: 1.73e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 3300, time: 22.917 days, wall time: 295.343 ms, max|u|: (2.90e-04, 2.38e-04), u★: 0.00 m s⁻¹, Q: 121.82 W m⁻², T₀: 16.59 ᵒC, extrema(T): (14.45, 17.36) ᵒC, S₀: 38.73 g/kg, e₀: 2.03e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 3400, time: 23.611 days, wall time: 295.760 ms, max|u|: (1.47e-04, 4.48e-04), u★: 0.00 m s⁻¹, Q: 113.12 W m⁻², T₀: 16.58 ᵒC, extrema(T): (14.45, 17.36) ᵒC, S₀: 38.73 g/kg, e₀: 1.69e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 3500, time: 24.306 days, wall time: 285.403 ms, max|u|: (7.13e-05, 2.37e-04), u★: 0.00 m s⁻¹, Q: 174.42 W m⁻², T₀: 16.57 ᵒC, extrema(T): (14.45, 17.36) ᵒC, S₀: 38.73 g/kg, e₀: 1.33e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 3600, time: 25 days, wall time: 284.990 ms, max|u|: (9.98e-04, 6.17e-04), u★: 0.00 m s⁻¹, Q: 171.39 W m⁻², T₀: 16.56 ᵒC, extrema(T): (14.45, 17.34) ᵒC, S₀: 38.74 g/kg, e₀: 2.38e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 3700, time: 25.694 days, wall time: 304.965 ms, max|u|: (2.40e-04, 9.44e-04), u★: 0.00 m s⁻¹, Q: 178.08 W m⁻², T₀: 16.55 ᵒC, extrema(T): (14.45, 17.34) ᵒC, S₀: 38.74 g/kg, e₀: 2.72e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 3800, time: 26.389 days, wall time: 305.672 ms, max|u|: (2.16e-04, 8.52e-04), u★: 0.00 m s⁻¹, Q: 236.47 W m⁻², T₀: 16.53 ᵒC, extrema(T): (14.45, 17.32) ᵒC, S₀: 38.74 g/kg, e₀: 7.81e-05 m² s⁻²
[ Info: Ocean Station Papa, iter: 3900, time: 27.083 days, wall time: 276.645 ms, max|u|: (4.06e-05, 5.38e-04), u★: 0.00 m s⁻¹, Q: 239.84 W m⁻², T₀: 16.52 ᵒC, extrema(T): (14.45, 17.30) ᵒC, S₀: 38.75 g/kg, e₀: 3.52e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 4000, time: 27.778 days, wall time: 292.019 ms, max|u|: (1.69e-04, 4.70e-04), u★: 0.00 m s⁻¹, Q: 228.07 W m⁻², T₀: 16.50 ᵒC, extrema(T): (14.45, 17.30) ᵒC, S₀: 38.75 g/kg, e₀: 3.26e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 4100, time: 28.472 days, wall time: 273.571 ms, max|u|: (3.71e-04, 4.19e-04), u★: 0.00 m s⁻¹, Q: 157.10 W m⁻², T₀: 16.49 ᵒC, extrema(T): (14.45, 17.28) ᵒC, S₀: 38.76 g/kg, e₀: 6.67e-05 m² s⁻²
[ Info: Ocean Station Papa, iter: 4200, time: 29.167 days, wall time: 278.779 ms, max|u|: (4.00e-04, 1.25e-04), u★: 0.00 m s⁻¹, Q: 130.91 W m⁻², T₀: 16.48 ᵒC, extrema(T): (14.45, 17.28) ᵒC, S₀: 38.76 g/kg, e₀: 2.16e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 4300, time: 29.861 days, wall time: 282.872 ms, max|u|: (4.67e-04, 1.05e-04), u★: 0.00 m s⁻¹, Q: 190.35 W m⁻², T₀: 16.47 ᵒC, extrema(T): (14.45, 17.25) ᵒC, S₀: 38.76 g/kg, e₀: 3.03e-04 m² s⁻²
[ Info: Simulation is stopping after running for 24.290 seconds.
[ Info: Simulation time 30 days equals or exceeds stop time 30 days.
[ Info: Drawing frame 1 of 241...
[ Info: Drawing frame 2 of 241...
[ Info: Drawing frame 3 of 241...
[ Info: Drawing frame 4 of 241...
[ Info: Drawing frame 5 of 241...
[ Info: Drawing frame 6 of 241...
[ Info: Drawing frame 7 of 241...
[ Info: Drawing frame 8 of 241...
[ Info: Drawing frame 9 of 241...
[ Info: Drawing frame 10 of 241...
[ Info: Drawing frame 11 of 241...
[ Info: Drawing frame 12 of 241...
[ Info: Drawing frame 13 of 241...
[ Info: Drawing frame 14 of 241...
[ Info: Drawing frame 15 of 241...
[ Info: Drawing frame 16 of 241...
[ Info: Drawing frame 17 of 241...
[ Info: Drawing frame 18 of 241...
[ Info: Drawing frame 19 of 241...
[ Info: Drawing frame 20 of 241...
[ Info: Drawing frame 21 of 241...
[ Info: Drawing frame 22 of 241...
[ Info: Drawing frame 23 of 241...
[ Info: Drawing frame 24 of 241...
[ Info: Drawing frame 25 of 241...
[ Info: Drawing frame 26 of 241...
[ Info: Drawing frame 27 of 241...
[ Info: Drawing frame 28 of 241...
[ Info: Drawing frame 29 of 241...
[ Info: Drawing frame 30 of 241...
[ Info: Drawing frame 31 of 241...
[ Info: Drawing frame 32 of 241...
[ Info: Drawing frame 33 of 241...
[ Info: Drawing frame 34 of 241...
[ Info: Drawing frame 35 of 241...
[ Info: Drawing frame 36 of 241...
[ Info: Drawing frame 37 of 241...
[ Info: Drawing frame 38 of 241...
[ Info: Drawing frame 39 of 241...
[ Info: Drawing frame 40 of 241...
[ Info: Drawing frame 41 of 241...
[ Info: Drawing frame 42 of 241...
[ Info: Drawing frame 43 of 241...
[ Info: Drawing frame 44 of 241...
[ Info: Drawing frame 45 of 241...
[ Info: Drawing frame 46 of 241...
[ Info: Drawing frame 47 of 241...
[ Info: Drawing frame 48 of 241...
[ Info: Drawing frame 49 of 241...
[ Info: Drawing frame 50 of 241...
[ Info: Drawing frame 51 of 241...
[ Info: Drawing frame 52 of 241...
[ Info: Drawing frame 53 of 241...
[ Info: Drawing frame 54 of 241...
[ Info: Drawing frame 55 of 241...
[ Info: Drawing frame 56 of 241...
[ Info: Drawing frame 57 of 241...
[ Info: Drawing frame 58 of 241...
[ Info: Drawing frame 59 of 241...
[ Info: Drawing frame 60 of 241...
[ Info: Drawing frame 61 of 241...
[ Info: Drawing frame 62 of 241...
[ Info: Drawing frame 63 of 241...
[ Info: Drawing frame 64 of 241...
[ Info: Drawing frame 65 of 241...
[ Info: Drawing frame 66 of 241...
[ Info: Drawing frame 67 of 241...
[ Info: Drawing frame 68 of 241...
[ Info: Drawing frame 69 of 241...
[ Info: Drawing frame 70 of 241...
[ Info: Drawing frame 71 of 241...
[ Info: Drawing frame 72 of 241...
[ Info: Drawing frame 73 of 241...
[ Info: Drawing frame 74 of 241...
[ Info: Drawing frame 75 of 241...
[ Info: Drawing frame 76 of 241...
[ Info: Drawing frame 77 of 241...
[ Info: Drawing frame 78 of 241...
[ Info: Drawing frame 79 of 241...
[ Info: Drawing frame 80 of 241...
[ Info: Drawing frame 81 of 241...
[ Info: Drawing frame 82 of 241...
[ Info: Drawing frame 83 of 241...
[ Info: Drawing frame 84 of 241...
[ Info: Drawing frame 85 of 241...
[ Info: Drawing frame 86 of 241...
[ Info: Drawing frame 87 of 241...
[ Info: Drawing frame 88 of 241...
[ Info: Drawing frame 89 of 241...
[ Info: Drawing frame 90 of 241...
[ Info: Drawing frame 91 of 241...
[ Info: Drawing frame 92 of 241...
[ Info: Drawing frame 93 of 241...
[ Info: Drawing frame 94 of 241...
[ Info: Drawing frame 95 of 241...
[ Info: Drawing frame 96 of 241...
[ Info: Drawing frame 97 of 241...
[ Info: Drawing frame 98 of 241...
[ Info: Drawing frame 99 of 241...
[ Info: Drawing frame 100 of 241...
[ Info: Drawing frame 101 of 241...
[ Info: Drawing frame 102 of 241...
[ Info: Drawing frame 103 of 241...
[ Info: Drawing frame 104 of 241...
[ Info: Drawing frame 105 of 241...
[ Info: Drawing frame 106 of 241...
[ Info: Drawing frame 107 of 241...
[ Info: Drawing frame 108 of 241...
[ Info: Drawing frame 109 of 241...
[ Info: Drawing frame 110 of 241...
[ Info: Drawing frame 111 of 241...
[ Info: Drawing frame 112 of 241...
[ Info: Drawing frame 113 of 241...
[ Info: Drawing frame 114 of 241...
[ Info: Drawing frame 115 of 241...
[ Info: Drawing frame 116 of 241...
[ Info: Drawing frame 117 of 241...
[ Info: Drawing frame 118 of 241...
[ Info: Drawing frame 119 of 241...
[ Info: Drawing frame 120 of 241...
[ Info: Drawing frame 121 of 241...
[ Info: Drawing frame 122 of 241...
[ Info: Drawing frame 123 of 241...
[ Info: Drawing frame 124 of 241...
[ Info: Drawing frame 125 of 241...
[ Info: Drawing frame 126 of 241...
[ Info: Drawing frame 127 of 241...
[ Info: Drawing frame 128 of 241...
[ Info: Drawing frame 129 of 241...
[ Info: Drawing frame 130 of 241...
[ Info: Drawing frame 131 of 241...
[ Info: Drawing frame 132 of 241...
[ Info: Drawing frame 133 of 241...
[ Info: Drawing frame 134 of 241...
[ Info: Drawing frame 135 of 241...
[ Info: Drawing frame 136 of 241...
[ Info: Drawing frame 137 of 241...
[ Info: Drawing frame 138 of 241...
[ Info: Drawing frame 139 of 241...
[ Info: Drawing frame 140 of 241...
[ Info: Drawing frame 141 of 241...
[ Info: Drawing frame 142 of 241...
[ Info: Drawing frame 143 of 241...
[ Info: Drawing frame 144 of 241...
[ Info: Drawing frame 145 of 241...
[ Info: Drawing frame 146 of 241...
[ Info: Drawing frame 147 of 241...
[ Info: Drawing frame 148 of 241...
[ Info: Drawing frame 149 of 241...
[ Info: Drawing frame 150 of 241...
[ Info: Drawing frame 151 of 241...
[ Info: Drawing frame 152 of 241...
[ Info: Drawing frame 153 of 241...
[ Info: Drawing frame 154 of 241...
[ Info: Drawing frame 155 of 241...
[ Info: Drawing frame 156 of 241...
[ Info: Drawing frame 157 of 241...
[ Info: Drawing frame 158 of 241...
[ Info: Drawing frame 159 of 241...
[ Info: Drawing frame 160 of 241...
[ Info: Drawing frame 161 of 241...
[ Info: Drawing frame 162 of 241...
[ Info: Drawing frame 163 of 241...
[ Info: Drawing frame 164 of 241...
[ Info: Drawing frame 165 of 241...
[ Info: Drawing frame 166 of 241...
[ Info: Drawing frame 167 of 241...
[ Info: Drawing frame 168 of 241...
[ Info: Drawing frame 169 of 241...
[ Info: Drawing frame 170 of 241...
[ Info: Drawing frame 171 of 241...
[ Info: Drawing frame 172 of 241...
[ Info: Drawing frame 173 of 241...
[ Info: Drawing frame 174 of 241...
[ Info: Drawing frame 175 of 241...
[ Info: Drawing frame 176 of 241...
[ Info: Drawing frame 177 of 241...
[ Info: Drawing frame 178 of 241...
[ Info: Drawing frame 179 of 241...
[ Info: Drawing frame 180 of 241...
[ Info: Drawing frame 181 of 241...
[ Info: Drawing frame 182 of 241...
[ Info: Drawing frame 183 of 241...
[ Info: Drawing frame 184 of 241...
[ Info: Drawing frame 185 of 241...
[ Info: Drawing frame 186 of 241...
[ Info: Drawing frame 187 of 241...
[ Info: Drawing frame 188 of 241...
[ Info: Drawing frame 189 of 241...
[ Info: Drawing frame 190 of 241...
[ Info: Drawing frame 191 of 241...
[ Info: Drawing frame 192 of 241...
[ Info: Drawing frame 193 of 241...
[ Info: Drawing frame 194 of 241...
[ Info: Drawing frame 195 of 241...
[ Info: Drawing frame 196 of 241...
[ Info: Drawing frame 197 of 241...
[ Info: Drawing frame 198 of 241...
[ Info: Drawing frame 199 of 241...
[ Info: Drawing frame 200 of 241...
[ Info: Drawing frame 201 of 241...
[ Info: Drawing frame 202 of 241...
[ Info: Drawing frame 203 of 241...
[ Info: Drawing frame 204 of 241...
[ Info: Drawing frame 205 of 241...
[ Info: Drawing frame 206 of 241...
[ Info: Drawing frame 207 of 241...
[ Info: Drawing frame 208 of 241...
[ Info: Drawing frame 209 of 241...
[ Info: Drawing frame 210 of 241...
[ Info: Drawing frame 211 of 241...
[ Info: Drawing frame 212 of 241...
[ Info: Drawing frame 213 of 241...
[ Info: Drawing frame 214 of 241...
[ Info: Drawing frame 215 of 241...
[ Info: Drawing frame 216 of 241...
[ Info: Drawing frame 217 of 241...
[ Info: Drawing frame 218 of 241...
[ Info: Drawing frame 219 of 241...
[ Info: Drawing frame 220 of 241...
[ Info: Drawing frame 221 of 241...
[ Info: Drawing frame 222 of 241...
[ Info: Drawing frame 223 of 241...
[ Info: Drawing frame 224 of 241...
[ Info: Drawing frame 225 of 241...
[ Info: Drawing frame 226 of 241...
[ Info: Drawing frame 227 of 241...
[ Info: Drawing frame 228 of 241...
[ Info: Drawing frame 229 of 241...
[ Info: Drawing frame 230 of 241...
[ Info: Drawing frame 231 of 241...
[ Info: Drawing frame 232 of 241...
[ Info: Drawing frame 233 of 241...
[ Info: Drawing frame 234 of 241...
[ Info: Drawing frame 235 of 241...
[ Info: Drawing frame 236 of 241...
[ Info: Drawing frame 237 of 241...
[ Info: Drawing frame 238 of 241...
[ Info: Drawing frame 239 of 241...
[ Info: Drawing frame 240 of 241...
[ Info: Drawing frame 241 of 241...
This page was generated using Literate.jl.