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 Oceananigans
using Oceananigans.Units
using Oceananigans.BuoyancyModels: buoyancy_frequency
using Oceananigans.Units: Time
using CairoMakie
using Printf
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
├── 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=ECCOMetadata(:temperature), S=ECCOMetadata(:salinity))
[ Info: Inpainting ClimaOcean.DataWrangling.ECCO.ECCO4Monthly temperature data from 1993-01-01T00:00:00...
[ Info: ... (1.224 minutes)
[ Info: Inpainting ClimaOcean.DataWrangling.ECCO.ECCO4Monthly salinity data from 1993-01-01T00:00:00...
[ Info: ... (1.242 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.
simulation_days = 31
snapshots_per_day = 8 # corresponding to JRA55's 3-hour frequency
last_time = simulation_days * snapshots_per_day
atmosphere = JRA55_prescribed_atmosphere(1:last_time;
longitude = λ★,
latitude = φ★,
backend = InMemory())
2×2×1×248 PrescribedAtmosphere{Float32} on Oceananigans.Grids.LatitudeLongitudeGrid:
├── times: 248-element StepRangeLen{Float64, Base.TwicePrecision{Float64}, Base.TwicePrecision{Float64}, Int64}
├── reference_height: 10.0
└── boundary_layer_height: 10.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 / days
fig = Figure(size=(800, 600))
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)
display(fig)
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.fluxes.total.ocean.momentum.u)
τy = first(sim.model.fluxes.total.ocean.momentum.v)
Q = first(sim.model.fluxes.total.ocean.heat)
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
end
simulation.callbacks[:progress] = Callback(progress, IterationInterval(100))
Callback of progress on IterationInterval(100)
Build flux outputs
τx = coupled_model.fluxes.total.ocean.momentum.u
τy = coupled_model.fluxes.total.ocean.momentum.v
JT = coupled_model.fluxes.total.ocean.tracers.T
Js = coupled_model.fluxes.total.ocean.tracers.S
E = coupled_model.fluxes.turbulent.fields.water_vapor
Qc = coupled_model.fluxes.turbulent.fields.sensible_heat
Qv = coupled_model.fluxes.turbulent.fields.latent_heat
ρₒ = coupled_model.fluxes.ocean_reference_density
cₚ = coupled_model.fluxes.ocean_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.4568, min=14.8668, mean=16.0725, 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.0664, min=38.5484, mean=38.9286, 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=4.09156e-5, min=0.0, mean=3.04304e-7)
Slice fields at the surface
outputs = merge(fields, fluxes)
filename = "single_column_omip_$(location_name)"
simulation.output_writers[:jld2] = JLD2OutputWriter(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
set_theme!(Theme(linewidth=3))
fig = Figure(size=(2400, 1800))
axτ = Axis(fig[1, 1:2], xlabel="Days since Oct 1 1992", ylabel="Wind stress (N m⁻²)")
axu = Axis(fig[2, 1:2], xlabel="Days since Oct 1 1992", ylabel="Velocities (m s⁻¹)")
axQ = Axis(fig[1, 3:4], xlabel="Days since Oct 1 1992", ylabel="Heat flux (W m⁻²)")
axT = Axis(fig[2, 3:4], xlabel="Days since Oct 1 1992", ylabel="Surface temperature (ᵒC)")
axF = Axis(fig[1, 5:6], xlabel="Days since Oct 1 1992", ylabel="Freshwater volume flux (m s⁻¹)")
axS = Axis(fig[2, 5:6], xlabel="Days since Oct 1 1992", ylabel="Surface salinity (g kg⁻¹)")
axuz = Axis(fig[3, 1], xlabel="Velocities (m s⁻¹)", ylabel="z (m)")
axTz = Axis(fig[3, 2], xlabel="Temperature (ᵒC)", ylabel="z (m)")
axSz = Axis(fig[3, 3], xlabel="Salinity (g kg⁻¹)", ylabel="z (m)")
axNz = Axis(fig[3, 4], xlabel="Buoyancy frequency (s⁻²)", ylabel="z (m)")
axκz = Axis(fig[3, 5], xlabel="Eddy diffusivity (m² s⁻¹)", ylabel="z (m)", xscale=log10)
axez = Axis(fig[3, 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)
slider = Slider(fig[4, 1:6], range=1:Nt, startvalue=1)
n = slider.value
times = (times .- times[1]) ./days
Nt = length(times)
tn = @lift times[$n]
colors = Makie.wong_colors()
ρₒ = coupled_model.fluxes.ocean_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)
Tmax = maximum(interior(T))
Tmin = minimum(interior(T))
xlims!(axTz, Tmin - 0.1, Tmax + 0.1)
Nmax = maximum(interior(N²))
Nmin = minimum(interior(N²))
xlims!(axNz, Nmin / 2, Nmin * 1.1)
emax = maximum(interior(e))
xlims!(axez, 8e-7, emax * 1.1)
xlims!(axκz, 1e-7, 10)
Smax = maximum(interior(S))
Smin = minimum(interior(S))
xlims!(axSz, Smin - 0.2, Smax + 0.2)
record(fig, "single_column_profiles.mp4", 1:2: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: 2.717 minutes, max|u|: (0.00e+00, 0.00e+00), u★: 0.00 m s⁻¹, Q: 219.42 W m⁻², T₀: 17.38 ᵒC, extrema(T): (14.87, 17.46) ᵒC, S₀: 38.55 g/kg, e₀: 0.00e+00 m² s⁻²
[ Info: ... simulation initialization complete (5.137 minutes)
[ Info: Executing initial time step...
┌ Warning: Simulation stopped during initialization.
└ @ Oceananigans.Simulations /central/scratch/esm/slurm-buildkite/climaocean-examples/1047/depot/default/packages/Oceananigans/dgGml/src/Simulations/run.jl:129
[ Info: ... initial time step complete (14.848 seconds).
[ Info: Ocean Station Papa, iter: 100, time: 16.667 hours, wall time: 5.675 minutes, max|u|: (1.53e-04, 7.66e-04), u★: 0.00 m s⁻¹, Q: 248.11 W m⁻², T₀: 17.34 ᵒC, extrema(T): (14.87, 17.46) ᵒC, S₀: 38.56 g/kg, e₀: 3.08e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 200, time: 1.389 days, wall time: 906.780 ms, max|u|: (3.11e-04, 1.91e-03), u★: 0.00 m s⁻¹, Q: 168.56 W m⁻², T₀: 17.30 ᵒC, extrema(T): (14.87, 17.46) ᵒC, S₀: 38.56 g/kg, e₀: 1.69e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 300, time: 2.083 days, wall time: 723.889 ms, max|u|: (1.72e-04, 9.31e-04), u★: 0.00 m s⁻¹, Q: 179.04 W m⁻², T₀: 17.27 ᵒC, extrema(T): (14.87, 17.46) ᵒC, S₀: 38.55 g/kg, e₀: 2.34e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 400, time: 2.778 days, wall time: 906.032 ms, max|u|: (1.31e-03, 1.29e-03), u★: 0.00 m s⁻¹, Q: 230.23 W m⁻², T₀: 17.25 ᵒC, extrema(T): (14.87, 17.46) ᵒC, S₀: 38.56 g/kg, e₀: 2.18e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 500, time: 3.472 days, wall time: 1.055 seconds, max|u|: (1.04e-03, 1.37e-03), u★: 0.00 m s⁻¹, Q: 242.34 W m⁻², T₀: 17.21 ᵒC, extrema(T): (14.87, 17.46) ᵒC, S₀: 38.56 g/kg, e₀: 2.91e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 600, time: 4.167 days, wall time: 1.053 seconds, max|u|: (7.69e-04, 1.44e-03), u★: 0.00 m s⁻¹, Q: 200.37 W m⁻², T₀: 17.17 ᵒC, extrema(T): (14.87, 17.46) ᵒC, S₀: 38.56 g/kg, e₀: 2.35e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 700, time: 4.861 days, wall time: 1.050 seconds, max|u|: (3.48e-04, 1.76e-03), u★: 0.00 m s⁻¹, Q: 163.72 W m⁻², T₀: 17.14 ᵒC, extrema(T): (14.87, 17.46) ᵒC, S₀: 38.56 g/kg, e₀: 2.09e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 800, time: 5.556 days, wall time: 1.124 seconds, max|u|: (8.22e-04, 1.54e-03), u★: 0.00 m s⁻¹, Q: 154.91 W m⁻², T₀: 17.12 ᵒC, extrema(T): (14.87, 17.45) ᵒC, S₀: 38.56 g/kg, e₀: 2.43e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 900, time: 6.250 days, wall time: 844.711 ms, max|u|: (2.20e-03, 1.45e-03), u★: 0.00 m s⁻¹, Q: 174.77 W m⁻², T₀: 17.09 ᵒC, extrema(T): (14.87, 17.36) ᵒC, S₀: 38.57 g/kg, e₀: 2.32e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 1000, time: 6.944 days, wall time: 1.426 seconds, max|u|: (1.90e-03, 2.58e-04), u★: 0.00 m s⁻¹, Q: 171.66 W m⁻², T₀: 17.06 ᵒC, extrema(T): (14.87, 17.36) ᵒC, S₀: 38.57 g/kg, e₀: 2.41e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 1100, time: 7.639 days, wall time: 981.665 ms, max|u|: (2.66e-03, 8.84e-04), u★: 0.00 m s⁻¹, Q: 185.65 W m⁻², T₀: 17.04 ᵒC, extrema(T): (14.87, 17.34) ᵒC, S₀: 38.57 g/kg, e₀: 2.37e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 1200, time: 8.333 days, wall time: 985.456 ms, max|u|: (1.89e-03, 5.69e-04), u★: 0.00 m s⁻¹, Q: 74.86 W m⁻², T₀: 17.01 ᵒC, extrema(T): (14.87, 17.27) ᵒC, S₀: 38.57 g/kg, e₀: 1.58e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 1300, time: 9.028 days, wall time: 1.922 seconds, max|u|: (1.91e-03, 1.16e-03), u★: 0.00 m s⁻¹, Q: 185.05 W m⁻², T₀: 16.98 ᵒC, extrema(T): (14.87, 17.27) ᵒC, S₀: 38.58 g/kg, e₀: 2.90e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 1400, time: 9.722 days, wall time: 1.434 seconds, max|u|: (1.28e-03, 2.04e-03), u★: 0.00 m s⁻¹, Q: 174.04 W m⁻², T₀: 16.96 ᵒC, extrema(T): (14.87, 17.27) ᵒC, S₀: 38.58 g/kg, e₀: 2.70e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 1500, time: 10.417 days, wall time: 1.838 seconds, max|u|: (3.37e-03, 4.44e-03), u★: 0.00 m s⁻¹, Q: 228.81 W m⁻², T₀: 16.92 ᵒC, extrema(T): (14.87, 17.18) ᵒC, S₀: 38.58 g/kg, e₀: 2.94e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 1600, time: 11.111 days, wall time: 1.774 seconds, max|u|: (1.70e-03, 3.15e-03), u★: 0.00 m s⁻¹, Q: 202.34 W m⁻², T₀: 16.89 ᵒC, extrema(T): (14.87, 17.17) ᵒC, S₀: 38.58 g/kg, e₀: 2.68e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 1700, time: 11.806 days, wall time: 2.518 seconds, max|u|: (2.02e-04, 3.08e-03), u★: 0.00 m s⁻¹, Q: 127.14 W m⁻², T₀: 16.87 ᵒC, extrema(T): (14.87, 17.17) ᵒC, S₀: 38.58 g/kg, e₀: 1.81e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 1800, time: 12.500 days, wall time: 1.764 seconds, max|u|: (1.60e-03, 3.17e-03), u★: 0.00 m s⁻¹, Q: 192.38 W m⁻², T₀: 16.85 ᵒC, extrema(T): (14.87, 17.16) ᵒC, S₀: 38.59 g/kg, e₀: 2.41e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 1900, time: 13.194 days, wall time: 2.316 seconds, max|u|: (2.37e-03, 4.02e-03), u★: 0.00 m s⁻¹, Q: 237.87 W m⁻², T₀: 16.81 ᵒC, extrema(T): (14.87, 17.08) ᵒC, S₀: 38.59 g/kg, e₀: 3.10e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 2000, time: 13.889 days, wall time: 2.843 seconds, max|u|: (1.56e-03, 3.36e-03), u★: 0.00 m s⁻¹, Q: 350.86 W m⁻², T₀: 16.77 ᵒC, extrema(T): (14.87, 17.08) ᵒC, S₀: 38.59 g/kg, e₀: 3.99e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 2100, time: 14.583 days, wall time: 1.974 seconds, max|u|: (2.03e-03, 2.71e-03), u★: 0.00 m s⁻¹, Q: 289.22 W m⁻², T₀: 16.73 ᵒC, extrema(T): (14.87, 17.02) ᵒC, S₀: 38.60 g/kg, e₀: 3.34e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 2200, time: 15.278 days, wall time: 2.386 seconds, max|u|: (1.98e-03, 1.57e-03), u★: 0.00 m s⁻¹, Q: 209.34 W m⁻², T₀: 16.68 ᵒC, extrema(T): (14.87, 16.98) ᵒC, S₀: 38.60 g/kg, e₀: 3.01e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 2300, time: 15.972 days, wall time: 2.378 seconds, max|u|: (2.30e-03, 8.67e-04), u★: 0.00 m s⁻¹, Q: 271.75 W m⁻², T₀: 16.64 ᵒC, extrema(T): (14.87, 16.98) ᵒC, S₀: 38.60 g/kg, e₀: 3.59e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 2400, time: 16.667 days, wall time: 2.136 seconds, max|u|: (1.90e-03, 1.45e-04), u★: 0.00 m s⁻¹, Q: 214.30 W m⁻², T₀: 16.61 ᵒC, extrema(T): (14.87, 16.91) ᵒC, S₀: 38.60 g/kg, e₀: 2.83e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 2500, time: 17.361 days, wall time: 2.041 seconds, max|u|: (1.88e-03, 1.27e-03), u★: 0.00 m s⁻¹, Q: 73.22 W m⁻², T₀: 16.58 ᵒC, extrema(T): (14.87, 16.90) ᵒC, S₀: 38.61 g/kg, e₀: 1.74e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 2600, time: 18.056 days, wall time: 2.426 seconds, max|u|: (1.15e-03, 1.60e-03), u★: 0.00 m s⁻¹, Q: 254.20 W m⁻², T₀: 16.55 ᵒC, extrema(T): (14.87, 16.90) ᵒC, S₀: 38.61 g/kg, e₀: 3.71e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 2700, time: 18.750 days, wall time: 1.729 seconds, max|u|: (2.41e-04, 1.42e-03), u★: 0.00 m s⁻¹, Q: 277.88 W m⁻², T₀: 16.53 ᵒC, extrema(T): (14.87, 16.90) ᵒC, S₀: 38.61 g/kg, e₀: 3.89e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 2800, time: 19.444 days, wall time: 1.739 seconds, max|u|: (2.87e-04, 1.26e-03), u★: 0.00 m s⁻¹, Q: 155.26 W m⁻², T₀: 16.49 ᵒC, extrema(T): (14.87, 16.82) ᵒC, S₀: 38.61 g/kg, e₀: 2.08e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 2900, time: 20.139 days, wall time: 1.612 seconds, max|u|: (8.22e-04, 1.10e-03), u★: 0.00 m s⁻¹, Q: 216.57 W m⁻², T₀: 16.46 ᵒC, extrema(T): (14.87, 16.82) ᵒC, S₀: 38.61 g/kg, e₀: 3.39e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 3000, time: 20.833 days, wall time: 1.510 seconds, max|u|: (1.65e-03, 1.07e-04), u★: 0.00 m s⁻¹, Q: 257.05 W m⁻², T₀: 16.43 ᵒC, extrema(T): (14.87, 16.82) ᵒC, S₀: 38.62 g/kg, e₀: 3.63e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 3100, time: 21.528 days, wall time: 1.590 seconds, max|u|: (1.51e-03, 5.89e-04), u★: 0.00 m s⁻¹, Q: 212.60 W m⁻², T₀: 16.41 ᵒC, extrema(T): (14.87, 16.82) ᵒC, S₀: 38.62 g/kg, e₀: 2.99e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 3200, time: 22.222 days, wall time: 2.114 seconds, max|u|: (9.40e-04, 1.02e-03), u★: 0.00 m s⁻¹, Q: 136.59 W m⁻², T₀: 16.38 ᵒC, extrema(T): (14.87, 16.74) ᵒC, S₀: 38.62 g/kg, e₀: 2.52e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 3300, time: 22.917 days, wall time: 1.398 seconds, max|u|: (5.27e-04, 1.82e-03), u★: 0.00 m s⁻¹, Q: 219.02 W m⁻², T₀: 16.36 ᵒC, extrema(T): (14.87, 16.74) ᵒC, S₀: 38.62 g/kg, e₀: 2.97e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 3400, time: 23.611 days, wall time: 1.445 seconds, max|u|: (5.86e-04, 2.15e-03), u★: 0.00 m s⁻¹, Q: 194.94 W m⁻², T₀: 16.33 ᵒC, extrema(T): (14.87, 16.74) ᵒC, S₀: 38.62 g/kg, e₀: 2.48e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 3500, time: 24.306 days, wall time: 1.473 seconds, max|u|: (1.05e-03, 1.18e-03), u★: 0.00 m s⁻¹, Q: 120.75 W m⁻², T₀: 16.30 ᵒC, extrema(T): (14.87, 16.71) ᵒC, S₀: 38.63 g/kg, e₀: 2.56e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 3600, time: 25 days, wall time: 1.454 seconds, max|u|: (2.96e-03, 2.71e-03), u★: 0.01 m s⁻¹, Q: 403.90 W m⁻², T₀: 16.26 ᵒC, extrema(T): (14.87, 16.66) ᵒC, S₀: 38.63 g/kg, e₀: 4.97e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 3700, time: 25.694 days, wall time: 1.643 seconds, max|u|: (9.05e-04, 2.52e-03), u★: 0.00 m s⁻¹, Q: 337.01 W m⁻², T₀: 16.22 ᵒC, extrema(T): (14.87, 16.66) ᵒC, S₀: 38.63 g/kg, e₀: 4.03e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 3800, time: 26.389 days, wall time: 1.541 seconds, max|u|: (5.68e-04, 1.79e-03), u★: 0.00 m s⁻¹, Q: 156.99 W m⁻², T₀: 16.18 ᵒC, extrema(T): (14.87, 16.58) ᵒC, S₀: 38.64 g/kg, e₀: 2.48e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 3900, time: 27.083 days, wall time: 1.188 seconds, max|u|: (2.04e-04, 1.57e-03), u★: 0.00 m s⁻¹, Q: 400.40 W m⁻², T₀: 16.13 ᵒC, extrema(T): (14.87, 16.58) ᵒC, S₀: 38.64 g/kg, e₀: 4.78e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 4000, time: 27.778 days, wall time: 1.862 seconds, max|u|: (1.55e-04, 1.99e-03), u★: 0.00 m s⁻¹, Q: 389.52 W m⁻², T₀: 16.09 ᵒC, extrema(T): (14.87, 16.53) ᵒC, S₀: 38.64 g/kg, e₀: 4.79e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 4100, time: 28.472 days, wall time: 1.374 seconds, max|u|: (1.03e-03, 1.57e-03), u★: 0.00 m s⁻¹, Q: 221.52 W m⁻², T₀: 16.05 ᵒC, extrema(T): (14.87, 16.53) ᵒC, S₀: 38.65 g/kg, e₀: 2.87e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 4200, time: 29.167 days, wall time: 1.723 seconds, max|u|: (1.23e-03, 7.02e-04), u★: 0.00 m s⁻¹, Q: 213.77 W m⁻², T₀: 16.02 ᵒC, extrema(T): (14.87, 16.52) ᵒC, S₀: 38.65 g/kg, e₀: 3.36e-04 m² s⁻²
[ Info: Ocean Station Papa, iter: 4300, time: 29.861 days, wall time: 1.372 seconds, max|u|: (1.99e-03, 2.26e-04), u★: 0.00 m s⁻¹, Q: 350.38 W m⁻², T₀: 15.99 ᵒC, extrema(T): (14.87, 16.49) ᵒC, S₀: 38.65 g/kg, e₀: 4.00e-04 m² s⁻²
[ Info: Simulation is stopping after running for 6.832 minutes.
[ Info: Simulation time 30 days equals or exceeds stop time 30 days.
[ Info: Drawing frame 1 of 241...
[ Info: Drawing frame 3 of 241...
[ Info: Drawing frame 5 of 241...
[ Info: Drawing frame 7 of 241...
[ Info: Drawing frame 9 of 241...
[ Info: Drawing frame 11 of 241...
[ Info: Drawing frame 13 of 241...
[ Info: Drawing frame 15 of 241...
[ Info: Drawing frame 17 of 241...
[ Info: Drawing frame 19 of 241...
[ Info: Drawing frame 21 of 241...
[ Info: Drawing frame 23 of 241...
[ Info: Drawing frame 25 of 241...
[ Info: Drawing frame 27 of 241...
[ Info: Drawing frame 29 of 241...
[ Info: Drawing frame 31 of 241...
[ Info: Drawing frame 33 of 241...
[ Info: Drawing frame 35 of 241...
[ Info: Drawing frame 37 of 241...
[ Info: Drawing frame 39 of 241...
[ Info: Drawing frame 41 of 241...
[ Info: Drawing frame 43 of 241...
[ Info: Drawing frame 45 of 241...
[ Info: Drawing frame 47 of 241...
[ Info: Drawing frame 49 of 241...
[ Info: Drawing frame 51 of 241...
[ Info: Drawing frame 53 of 241...
[ Info: Drawing frame 55 of 241...
[ Info: Drawing frame 57 of 241...
[ Info: Drawing frame 59 of 241...
[ Info: Drawing frame 61 of 241...
[ Info: Drawing frame 63 of 241...
[ Info: Drawing frame 65 of 241...
[ Info: Drawing frame 67 of 241...
[ Info: Drawing frame 69 of 241...
[ Info: Drawing frame 71 of 241...
[ Info: Drawing frame 73 of 241...
[ Info: Drawing frame 75 of 241...
[ Info: Drawing frame 77 of 241...
[ Info: Drawing frame 79 of 241...
[ Info: Drawing frame 81 of 241...
[ Info: Drawing frame 83 of 241...
[ Info: Drawing frame 85 of 241...
[ Info: Drawing frame 87 of 241...
[ Info: Drawing frame 89 of 241...
[ Info: Drawing frame 91 of 241...
[ Info: Drawing frame 93 of 241...
[ Info: Drawing frame 95 of 241...
[ Info: Drawing frame 97 of 241...
[ Info: Drawing frame 99 of 241...
[ Info: Drawing frame 101 of 241...
[ Info: Drawing frame 103 of 241...
[ Info: Drawing frame 105 of 241...
[ Info: Drawing frame 107 of 241...
[ Info: Drawing frame 109 of 241...
[ Info: Drawing frame 111 of 241...
[ Info: Drawing frame 113 of 241...
[ Info: Drawing frame 115 of 241...
[ Info: Drawing frame 117 of 241...
[ Info: Drawing frame 119 of 241...
[ Info: Drawing frame 121 of 241...
[ Info: Drawing frame 123 of 241...
[ Info: Drawing frame 125 of 241...
[ Info: Drawing frame 127 of 241...
[ Info: Drawing frame 129 of 241...
[ Info: Drawing frame 131 of 241...
[ Info: Drawing frame 133 of 241...
[ Info: Drawing frame 135 of 241...
[ Info: Drawing frame 137 of 241...
[ Info: Drawing frame 139 of 241...
[ Info: Drawing frame 141 of 241...
[ Info: Drawing frame 143 of 241...
[ Info: Drawing frame 145 of 241...
[ Info: Drawing frame 147 of 241...
[ Info: Drawing frame 149 of 241...
[ Info: Drawing frame 151 of 241...
[ Info: Drawing frame 153 of 241...
[ Info: Drawing frame 155 of 241...
[ Info: Drawing frame 157 of 241...
[ Info: Drawing frame 159 of 241...
[ Info: Drawing frame 161 of 241...
[ Info: Drawing frame 163 of 241...
[ Info: Drawing frame 165 of 241...
[ Info: Drawing frame 167 of 241...
[ Info: Drawing frame 169 of 241...
[ Info: Drawing frame 171 of 241...
[ Info: Drawing frame 173 of 241...
[ Info: Drawing frame 175 of 241...
[ Info: Drawing frame 177 of 241...
[ Info: Drawing frame 179 of 241...
[ Info: Drawing frame 181 of 241...
[ Info: Drawing frame 183 of 241...
[ Info: Drawing frame 185 of 241...
[ Info: Drawing frame 187 of 241...
[ Info: Drawing frame 189 of 241...
[ Info: Drawing frame 191 of 241...
[ Info: Drawing frame 193 of 241...
[ Info: Drawing frame 195 of 241...
[ Info: Drawing frame 197 of 241...
[ Info: Drawing frame 199 of 241...
[ Info: Drawing frame 201 of 241...
[ Info: Drawing frame 203 of 241...
[ Info: Drawing frame 205 of 241...
[ Info: Drawing frame 207 of 241...
[ Info: Drawing frame 209 of 241...
[ Info: Drawing frame 211 of 241...
[ Info: Drawing frame 213 of 241...
[ Info: Drawing frame 215 of 241...
[ Info: Drawing frame 217 of 241...
[ Info: Drawing frame 219 of 241...
[ Info: Drawing frame 221 of 241...
[ Info: Drawing frame 223 of 241...
[ Info: Drawing frame 225 of 241...
[ Info: Drawing frame 227 of 241...
[ Info: Drawing frame 229 of 241...
[ Info: Drawing frame 231 of 241...
[ Info: Drawing frame 233 of 241...
[ Info: Drawing frame 235 of 241...
[ Info: Drawing frame 237 of 241...
[ Info: Drawing frame 239 of 241...
[ Info: Drawing frame 241 of 241...
This page was generated using Literate.jl.