Wind- and convection-driven mixing in an ocean surface boundary layer

This example simulates mixing by three-dimensional turbulence in an ocean surface boundary layer driven by atmospheric winds and convection. It demonstrates:

  • How to set-up a grid with varying spacing in the vertical direction
  • How to use the SeawaterBuoyancy model for buoyancy with a linear equation of state.
  • How to use a turbulence closure for large eddy simulation.
  • How to use a function to impose a boundary condition.

Install dependencies

First let's make sure we have all required packages installed.

using Pkg
pkg"add Oceananigans, JLD2, Plots"

We start by importing all of the packages and functions that we'll need for this example.

using Random
using Printf
using Plots
using JLD2

using Oceananigans
using Oceananigans.Units: minute, minutes, hour

The grid

We use 32³ grid points with 2 m grid spacing in the horizontal and varying spacing in the vertical, with higher resolution closer to the surface,

σ = 1.1 # stretching factor
Nz = 24
Lz = 32

hyperbolically_spaced_faces(k) = - Lz * (1 - tanh(σ * (k - 1) / Nz) / tanh(σ))

grid = VerticallyStretchedRectilinearGrid(size = (32, 32, Nz),
                                          x = (0, 64),
                                          y = (0, 64),
                                          z_faces = hyperbolically_spaced_faces)
VerticallyStretchedRectilinearGrid{Float64, Periodic, Periodic, Bounded}
                   domain: x ∈ [0.0, 64.0], y ∈ [0.0, 64.0], z ∈ [-32.0, -0.0]
                 topology: (Periodic, Periodic, Bounded)
  resolution (Nx, Ny, Nz): (32, 32, 24)
   halo size (Hx, Hy, Hz): (1, 1, 1)
grid spacing (Δx, Δy, Δz): (2.0, 2.0, [min=0.6826950100338962, max=1.8309085743885056])

We plot vertical spacing versus depth to inspect the prescribed grid stretching:

plot(grid.Δzᵃᵃᶜ[1:Nz], grid.zᵃᵃᶜ[1:Nz],
      marker = :circle,
      ylabel = "Depth (m)",
      xlabel = "Vertical spacing (m)",
      legend = nothing)

Buoyancy that depends on temperature and salinity

We use the SeawaterBuoyancy model with a linear equation of state,

buoyancy = SeawaterBuoyancy(equation_of_state=LinearEquationOfState(α=2e-4, β=8e-4))
SeawaterBuoyancy{Float64}: g = 9.80665
└── equation of state: LinearEquationOfState{Float64}: α = 2.00e-04, β = 8.00e-04

where $α$ and $β$ are the thermal expansion and haline contraction coefficients for temperature and salinity.

Boundary conditions

We calculate the surface temperature flux associated with surface heating of 200 W m⁻², reference density ρₒ, and heat capacity cᴾ,

Qʰ = 200  # W m⁻², surface _heat_ flux
ρₒ = 1026 # kg m⁻³, average density at the surface of the world ocean
cᴾ = 3991 # J K⁻¹ kg⁻¹, typical heat capacity for seawater

Qᵀ = Qʰ / (ρₒ * cᴾ) # K m s⁻¹, surface _temperature_ flux
4.884283985946938e-5

Finally, we impose a temperature gradient dTdz both initially and at the bottom of the domain, culminating in the boundary conditions on temperature,

dTdz = 0.01 # K m⁻¹

T_bcs = FieldBoundaryConditions(top = FluxBoundaryCondition(Qᵀ),
                                bottom = GradientBoundaryCondition(dTdz))
Oceananigans.FieldBoundaryConditions, with boundary conditions
├── west: Oceananigans.BoundaryConditions.DefaultPrognosticFieldBoundaryCondition
├── east: Oceananigans.BoundaryConditions.DefaultPrognosticFieldBoundaryCondition
├── south: Oceananigans.BoundaryConditions.DefaultPrognosticFieldBoundaryCondition
├── north: Oceananigans.BoundaryConditions.DefaultPrognosticFieldBoundaryCondition
├── bottom: BoundaryCondition{Gradient, Float64}
├── top: BoundaryCondition{Flux, Float64}
└── immersed: BoundaryCondition{Flux, Nothing}

Note that a positive temperature flux at the surface of the ocean implies cooling. This is because a positive temperature flux implies that temperature is fluxed upwards, out of the ocean.

For the velocity field, we imagine a wind blowing over the ocean surface with an average velocity at 10 meters u₁₀, and use a drag coefficient cᴰ to estimate the kinematic stress (that is, stress divided by density) exerted by the wind on the ocean:

u₁₀ = 10     # m s⁻¹, average wind velocity 10 meters above the ocean
cᴰ = 2.5e-3  # dimensionless drag coefficient
ρₐ = 1.225   # kg m⁻³, average density of air at sea-level

Qᵘ = - ρₐ / ρₒ * cᴰ * u₁₀ * abs(u₁₀) # m² s⁻²
-0.0002984892787524367

The boundary conditions on u are thus

u_bcs = FieldBoundaryConditions(top = FluxBoundaryCondition(Qᵘ))
Oceananigans.FieldBoundaryConditions, with boundary conditions
├── west: Oceananigans.BoundaryConditions.DefaultPrognosticFieldBoundaryCondition
├── east: Oceananigans.BoundaryConditions.DefaultPrognosticFieldBoundaryCondition
├── south: Oceananigans.BoundaryConditions.DefaultPrognosticFieldBoundaryCondition
├── north: Oceananigans.BoundaryConditions.DefaultPrognosticFieldBoundaryCondition
├── bottom: Oceananigans.BoundaryConditions.DefaultPrognosticFieldBoundaryCondition
├── top: BoundaryCondition{Flux, Float64}
└── immersed: BoundaryCondition{Flux, Nothing}

For salinity, S, we impose an evaporative flux of the form

@inline Qˢ(x, y, t, S, evaporation_rate) = - evaporation_rate * S # [salinity unit] m s⁻¹

where S is salinity. We use an evporation rate of 1 millimeter per hour,

evaporation_rate = 1e-3 / hour # m s⁻¹
2.7777777777777776e-7

We build the Flux evaporation BoundaryCondition with the function , indicating that depends on salinity S and passing the parameter evaporation_rate,

evaporation_bc = FluxBoundaryCondition(Qˢ, field_dependencies=:S, parameters=evaporation_rate)
BoundaryCondition: classification=Flux, condition=Qˢ(x, y, t, S, evaporation_rate) in Main at ocean_wind_mixing_and_convection.md:126

The full salinity boundary conditions are

S_bcs = FieldBoundaryConditions(top=evaporation_bc)
Oceananigans.FieldBoundaryConditions, with boundary conditions
├── west: Oceananigans.BoundaryConditions.DefaultPrognosticFieldBoundaryCondition
├── east: Oceananigans.BoundaryConditions.DefaultPrognosticFieldBoundaryCondition
├── south: Oceananigans.BoundaryConditions.DefaultPrognosticFieldBoundaryCondition
├── north: Oceananigans.BoundaryConditions.DefaultPrognosticFieldBoundaryCondition
├── bottom: Oceananigans.BoundaryConditions.DefaultPrognosticFieldBoundaryCondition
├── top: BoundaryCondition{Flux, Oceananigans.BoundaryConditions.ContinuousBoundaryFunction{Nothing, Nothing, Nothing, Nothing, typeof(Main.Qˢ), Float64, Tuple{Symbol}, Nothing, Nothing}}
└── immersed: BoundaryCondition{Flux, Nothing}

Model instantiation

We fill in the final details of the model here: upwind-biased 5th-order advection for momentum and tracers, 3rd-order Runge-Kutta time-stepping, Coriolis forces, and the AnisotropicMinimumDissipation closure for large eddy simulation to model the effect of turbulent motions at scales smaller than the grid scale that we cannot explicitly resolve.

model = NonhydrostaticModel(architecture = CPU(),
                            advection = UpwindBiasedFifthOrder(),
                            timestepper = :RungeKutta3,
                            grid = grid,
                            tracers = (:T, :S),
                            coriolis = FPlane(f=1e-4),
                            buoyancy = buoyancy,
                            closure = AnisotropicMinimumDissipation(),
                            boundary_conditions = (u=u_bcs, T=T_bcs, S=S_bcs))
NonhydrostaticModel{CPU, Float64}(time = 0 seconds, iteration = 0) 
├── grid: VerticallyStretchedRectilinearGrid{Float64, Periodic, Periodic, Bounded}(Nx=32, Ny=32, Nz=24)
├── tracers: (:T, :S)
├── closure: AnisotropicMinimumDissipation{Oceananigans.TurbulenceClosures.ExplicitTimeDiscretization, Float64, NamedTuple{(:T, :S), Tuple{Float64, Float64}}, Float64, NamedTuple{(:T, :S), Tuple{Float64, Float64}}, Nothing}
├── buoyancy: SeawaterBuoyancy{Float64, LinearEquationOfState{Float64}, Nothing, Nothing}
└── coriolis: FPlane{Float64}

Notes:

  • To use the Smagorinsky-Lilly turbulence closure (with a constant model coefficient) rather than AnisotropicMinimumDissipation, use closure = SmagorinskyLilly() in the model constructor.

  • To change the architecture to GPU, replace architecture = CPU() with architecture = GPU().

Initial conditions

Our initial condition for temperature consists of a linear stratification superposed with random noise damped at the walls, while our initial condition for velocity consists only of random noise.

# Random noise damped at top and bottom
Ξ(z) = randn() * z / model.grid.Lz * (1 + z / model.grid.Lz) # noise

# Temperature initial condition: a stable density gradient with random noise superposed.
Tᵢ(x, y, z) = 20 + dTdz * z + dTdz * model.grid.Lz * 1e-6 * Ξ(z)

# Velocity initial condition: random noise scaled by the friction velocity.
uᵢ(x, y, z) = sqrt(abs(Qᵘ)) * 1e-3 * Ξ(z)

# `set!` the `model` fields using functions or constants:
set!(model, u=uᵢ, w=uᵢ, T=Tᵢ, S=35)

Setting up a simulation

We first build a TimeStepWizard to ensure stable time-stepping with a Courant-Freidrichs-Lewy (CFL) number of 1.0.

wizard = TimeStepWizard(cfl=1.0, Δt=10.0, max_change=1.1, max_Δt=1minute)
TimeStepWizard{Float64, typeof(Oceananigans.Utils.cell_advection_timescale), typeof(Oceananigans.Simulations.infinite_diffusion_timescale)}(1.0, Inf, 1.1, 0.5, 60.0, 0.0, 10.0, Oceananigans.Utils.cell_advection_timescale, Oceananigans.Simulations.infinite_diffusion_timescale)

Nice progress messaging is helpful:

start_time = time_ns() # so we can print the total elapsed wall time

# Print a progress message
progress_message(sim) =
    @printf("i: %04d, t: %s, Δt: %s, wmax = %.1e ms⁻¹, wall time: %s\n",
            sim.model.clock.iteration, prettytime(model.clock.time),
            prettytime(wizard.Δt), maximum(abs, sim.model.velocities.w),
            prettytime((time_ns() - start_time) * 1e-9))
progress_message (generic function with 1 method)

We then set up the simulation:

simulation = Simulation(model, Δt=wizard, stop_time=40minutes, iteration_interval=10,
                        progress=progress_message)
Simulation{typename(NonhydrostaticModel){typename(CPU), Float64}}
├── Model clock: time = 0 seconds, iteration = 0 
├── Next time step (TimeStepWizard{Float64, typeof(Oceananigans.Utils.cell_advection_timescale), typeof(Oceananigans.Simulations.infinite_diffusion_timescale)}): 10 seconds 
├── Iteration interval: 10
├── Stop criteria: Any[Oceananigans.Simulations.iteration_limit_exceeded, Oceananigans.Simulations.stop_time_exceeded, Oceananigans.Simulations.wall_time_limit_exceeded]
├── Run time: 0 seconds, wall time limit: Inf
├── Stop time: 40 minutes, stop iteration: Inf
├── Diagnostics: typename(OrderedCollections.OrderedDict) with 1 entry:
│   └── nan_checker => typename(NaNChecker)
└── Output writers: typename(OrderedCollections.OrderedDict) with no entries

Output

We use the JLD2OutputWriter to save $x, z$ slices of the velocity fields, tracer fields, and eddy diffusivities. The prefix keyword argument to JLD2OutputWriter indicates that output will be saved in ocean_wind_mixing_and_convection.jld2.

# Create a NamedTuple with eddy viscosity
eddy_viscosity = (νₑ = model.diffusivity_fields.νₑ,)

simulation.output_writers[:slices] =
    JLD2OutputWriter(model, merge(model.velocities, model.tracers, eddy_viscosity),
                           prefix = "ocean_wind_mixing_and_convection",
                     field_slicer = FieldSlicer(j=Int(grid.Ny/2)),
                         schedule = TimeInterval(1minute),
                            force = true)
JLD2OutputWriter scheduled on TimeInterval(1 minute):
├── filepath: ./ocean_wind_mixing_and_convection.jld2
├── 6 outputs: (:u, :v, :w, :T, :S, :νₑ)
├── field slicer: FieldSlicer(:, 16, :, with_halos=false)
├── array type: Array{Float32}
├── including: [:grid, :coriolis, :buoyancy, :closure]
└── max filesize: Inf YiB

We're ready:

run!(simulation)
[ Info: Updating model auxiliary state before the first time step...
[ Info:     ... updated in 7.716 ms.
[ Info: Executing first time step...
i: 0010, t: 1.667 minutes, Δt: 10 seconds, wmax = 1.1e-05 ms⁻¹, wall time: 1.434 minutes
i: 0020, t: 3 minutes, Δt: 11 seconds, wmax = 9.9e-06 ms⁻¹, wall time: 1.449 minutes
i: 0030, t: 5 minutes, Δt: 12.100 seconds, wmax = 7.3e-06 ms⁻¹, wall time: 1.460 minutes
i: 0040, t: 7 minutes, Δt: 13.310 seconds, wmax = 6.6e-06 ms⁻¹, wall time: 1.472 minutes
i: 0050, t: 8.726 minutes, Δt: 10.894 seconds, wmax = 3.6e-05 ms⁻¹, wall time: 1.483 minutes
i: 0060, t: 10.146 minutes, Δt: 8.738 seconds, wmax = 2.6e-04 ms⁻¹, wall time: 1.493 minutes
i: 0070, t: 11.375 minutes, Δt: 7.504 seconds, wmax = 1.4e-03 ms⁻¹, wall time: 1.503 minutes
i: 0080, t: 12.442 minutes, Δt: 6.637 seconds, wmax = 7.0e-03 ms⁻¹, wall time: 1.512 minutes
i: 0090, t: 13.396 minutes, Δt: 5.947 seconds, wmax = 1.8e-02 ms⁻¹, wall time: 1.524 minutes
i: 0100, t: 14.269 minutes, Δt: 5.387 seconds, wmax = 3.2e-02 ms⁻¹, wall time: 1.535 minutes
i: 0110, t: 15.086 minutes, Δt: 5.173 seconds, wmax = 4.7e-02 ms⁻¹, wall time: 1.546 minutes
i: 0120, t: 15.920 minutes, Δt: 5.002 seconds, wmax = 4.7e-02 ms⁻¹, wall time: 1.557 minutes
i: 0130, t: 16.825 minutes, Δt: 5.502 seconds, wmax = 5.1e-02 ms⁻¹, wall time: 1.567 minutes
i: 0140, t: 17.790 minutes, Δt: 5.926 seconds, wmax = 5.6e-02 ms⁻¹, wall time: 1.578 minutes
i: 0150, t: 18.848 minutes, Δt: 6.363 seconds, wmax = 5.4e-02 ms⁻¹, wall time: 1.589 minutes
i: 0160, t: 19.781 minutes, Δt: 6.691 seconds, wmax = 5.5e-02 ms⁻¹, wall time: 1.600 minutes
i: 0170, t: 20.811 minutes, Δt: 6.950 seconds, wmax = 6.0e-02 ms⁻¹, wall time: 1.610 minutes
i: 0180, t: 21.925 minutes, Δt: 6.940 seconds, wmax = 5.8e-02 ms⁻¹, wall time: 1.620 minutes
i: 0190, t: 23 minutes, Δt: 6.919 seconds, wmax = 6.4e-02 ms⁻¹, wall time: 1.631 minutes
i: 0200, t: 24.119 minutes, Δt: 7.142 seconds, wmax = 6.4e-02 ms⁻¹, wall time: 1.641 minutes
i: 0210, t: 25.225 minutes, Δt: 6.759 seconds, wmax = 5.5e-02 ms⁻¹, wall time: 1.651 minutes
i: 0220, t: 26.213 minutes, Δt: 6.403 seconds, wmax = 4.8e-02 ms⁻¹, wall time: 1.663 minutes
i: 0230, t: 27.224 minutes, Δt: 6.731 seconds, wmax = 4.8e-02 ms⁻¹, wall time: 1.673 minutes
i: 0240, t: 28.351 minutes, Δt: 7.024 seconds, wmax = 4.7e-02 ms⁻¹, wall time: 1.685 minutes
i: 0250, t: 29.459 minutes, Δt: 6.892 seconds, wmax = 4.9e-02 ms⁻¹, wall time: 1.697 minutes
i: 0260, t: 30.548 minutes, Δt: 6.574 seconds, wmax = 5.3e-02 ms⁻¹, wall time: 1.709 minutes
i: 0270, t: 31.685 minutes, Δt: 6.852 seconds, wmax = 4.7e-02 ms⁻¹, wall time: 1.718 minutes
i: 0280, t: 32.785 minutes, Δt: 6.725 seconds, wmax = 5.2e-02 ms⁻¹, wall time: 1.730 minutes
i: 0290, t: 33.899 minutes, Δt: 6.743 seconds, wmax = 5.0e-02 ms⁻¹, wall time: 1.743 minutes
i: 0300, t: 34.944 minutes, Δt: 6.294 seconds, wmax = 5.2e-02 ms⁻¹, wall time: 1.753 minutes
i: 0310, t: 35.986 minutes, Δt: 6.574 seconds, wmax = 4.7e-02 ms⁻¹, wall time: 1.762 minutes
i: 0320, t: 37 minutes, Δt: 6.693 seconds, wmax = 4.6e-02 ms⁻¹, wall time: 1.774 minutes
i: 0330, t: 38.113 minutes, Δt: 6.763 seconds, wmax = 3.7e-02 ms⁻¹, wall time: 1.783 minutes
i: 0340, t: 39.104 minutes, Δt: 6.245 seconds, wmax = 4.2e-02 ms⁻¹, wall time: 1.793 minutes
i: 0349, t: 40 minutes, Δt: 6.443 seconds, wmax = 4.4e-02 ms⁻¹, wall time: 1.804 minutes
[ Info: Simulation is stopping. Model time 40 minutes has hit or exceeded simulation stop time 40 minutes.

Turbulence visualization

We animate the data saved in ocean_wind_mixing_and_convection.jld2. We prepare for animating the flow by creating coordinate arrays, opening the file, building a vector of the iterations that we saved data at, and defining functions for computing colorbar limits:

# Coordinate arrays
xw, yw, zw = nodes(model.velocities.w)
xT, yT, zT = nodes(model.tracers.T)

# Open the file with our data
file = jldopen(simulation.output_writers[:slices].filepath)

# Extract a vector of iterations
iterations = parse.(Int, keys(file["timeseries/t"]))

""" Returns colorbar levels equispaced between `(-clim, clim)` and encompassing the extrema of `c`. """
function divergent_levels(c, clim, nlevels=21)
    cmax = maximum(abs, c)
    levels = clim > cmax ? range(-clim, stop=clim, length=nlevels) : range(-cmax, stop=cmax, length=nlevels)
    return (levels[1], levels[end]), levels
end

""" Returns colorbar levels equispaced between `clims` and encompassing the extrema of `c`."""
function sequential_levels(c, clims, nlevels=20)
    levels = range(clims[1], stop=clims[2], length=nlevels)
    cmin, cmax = minimum(c), maximum(c)
    cmin < clims[1] && (levels = vcat([cmin], levels))
    cmax > clims[2] && (levels = vcat(levels, [cmax]))
    return clims, levels
end

We start the animation at t = 10minutes since things are pretty boring till then:

times = [file["timeseries/t/$iter"] for iter in iterations]
intro = searchsortedfirst(times, 10minutes)

anim = @animate for (i, iter) in enumerate(iterations[intro:end])

    @info "Drawing frame $i from iteration $iter..."

    t = file["timeseries/t/$iter"]
    w = file["timeseries/w/$iter"][:, 1, :]
    T = file["timeseries/T/$iter"][:, 1, :]
    S = file["timeseries/S/$iter"][:, 1, :]
    νₑ = file["timeseries/νₑ/$iter"][:, 1, :]

    wlims, wlevels = divergent_levels(w, 2e-2)
    Tlims, Tlevels = sequential_levels(T, (19.7, 19.99))
    Slims, Slevels = sequential_levels(S, (35, 35.005))
    νlims, νlevels = sequential_levels(νₑ, (1e-6, 5e-3))

    kwargs = (linewidth=0, xlabel="x (m)", ylabel="z (m)", aspectratio=1,
              xlims=(0, grid.Lx), ylims=(-grid.Lz, 0))

    w_plot = contourf(xw, zw, w'; color=:balance, clims=wlims, levels=wlevels, kwargs...)
    T_plot = contourf(xT, zT, T'; color=:thermal, clims=Tlims, levels=Tlevels, kwargs...)
    S_plot = contourf(xT, zT, S'; color=:haline,  clims=Slims, levels=Slevels, kwargs...)

    # We use a heatmap for the eddy viscosity to observe how it varies on the grid scale.
    ν_plot = heatmap(xT, zT, νₑ'; color=:thermal, clims=νlims, levels=νlevels, kwargs...)

    w_title = @sprintf("vertical velocity (m s⁻¹), t = %s", prettytime(t))
    T_title = "temperature (ᵒC)"
    S_title = "salinity (g kg⁻¹)"
    ν_title = "eddy viscosity (m² s⁻¹)"

    # Arrange the plots side-by-side.
    plot(w_plot, T_plot, S_plot, ν_plot, layout=(2, 2), size=(1200, 600),
         title=[w_title T_title S_title ν_title])

    iter == iterations[end] && close(file)
end

This page was generated using Literate.jl.