Langmuir turbulence example

This example implements the Langmuir turbulence simulation reported in section 4 of

McWilliams, J. C. et al., "Langmuir Turbulence in the ocean," Journal of Fluid Mechanics (1997).

This example demonstrates

  • How to run large eddy simulations with surface wave effects via the Craik-Leibovich approximation.

  • How to specify time- and horizontally-averaged output.

Install dependencies

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

using Pkg
pkg"add Oceananigans, JLD2, Plots"
using Oceananigans

Model set-up

To build the model, we specify the grid, Stokes drift, boundary conditions, and Coriolis parameter.

Domain and numerical grid specification

We create a grid with modest resolution. The grid extent is similar, but not exactly the same as that in McWilliams et al. (1997).

grid = RegularCartesianGrid(size=(32, 32, 48), extent=(128, 128, 96))
RegularCartesianGrid{Float64, Periodic, Periodic, Bounded}
                   domain: x ∈ [0.0, 128.0], y ∈ [0.0, 128.0], z ∈ [-96.0, 0.0]
                 topology: (Periodic, Periodic, Bounded)
  resolution (Nx, Ny, Nz): (32, 32, 48)
   halo size (Hx, Hy, Hz): (1, 1, 1)
grid spacing (Δx, Δy, Δz): (4.0, 4.0, 2.0)

The Stokes Drift profile

The surface wave Stokes drift profile prescribed in McWilliams et al. (1997) corresponds to a 'monochromatic' (that is, single-frequency) wave field.

A monochromatic wave field is characterized by its wavelength and amplitude (half the distance from wave crest to wave trough), which determine the wave frequency and the vertical scale of the Stokes drift profile.

using Oceananigans.Buoyancy: g_Earth

 amplitude = 0.8 # m
wavelength = 60 # m
wavenumber = 2π / wavelength # m⁻¹
 frequency = sqrt(g_Earth * wavenumber) # s⁻¹

# The vertical scale over which the Stokes drift of a monochromatic surface wave
# decays away from the surface is `1/2wavenumber`, or
const vertical_scale = wavelength / 4π

# Stokes drift velocity at the surface
const Uˢ = amplitude^2 * wavenumber * frequency # m s⁻¹
0.06791774197745354

The const declarations ensure that Stokes drift functions compile on the GPU. To run this example on the GPU, write architecture = GPU() in the constructor for IncompressibleModel below.

The Stokes drift profile is

uˢ(z) = Uˢ * exp(z / vertical_scale)
uˢ (generic function with 1 method)

which we'll need for the initial condition.

The Craik-Leibovich equations in Oceananigans

Oceananigans implements the Craik-Leibovich approximation for surface wave effects using the Lagrangian-mean velocity field as its prognostic momentum variable. In other words, model.velocities.u is the Lagrangian-mean $x$-velocity beneath surface waves. This differs from models that use the Eulerian-mean velocity field as a prognostic variable, but has the advantage that $u$ accounts for the total advection of tracers and momentum, and that $u = v = w = 0$ is a steady solution even when Coriolis forces are present. See the physics documentation for more information.

The vertical derivative of the Stokes drift is

∂z_uˢ(z, t) = 1 / vertical_scale * Uˢ * exp(z / vertical_scale)
∂z_uˢ (generic function with 1 method)

Finally, we note that the time-derivative of the Stokes drift must be provided if the Stokes drift changes in time. In this example, the Stokes drift is constant and thus the time-derivative of the Stokes drift is 0.

Boundary conditions

At the surface at $z=0$, McWilliams et al. (1997) impose a wind stress on $u$,

using Oceananigans.BoundaryConditions

Qᵘ = -3.72e-5 # m² s⁻², surface kinematic momentum flux

u_boundary_conditions = UVelocityBoundaryConditions(grid, top = BoundaryCondition(Flux, Qᵘ))
Oceananigans.FieldBoundaryConditions (NamedTuple{(:x, :y, :z)}), with boundary conditions
├── x: CoordinateBoundaryConditions{BoundaryCondition{Oceananigans.BoundaryConditions.Periodic,Nothing},BoundaryCondition{Oceananigans.BoundaryConditions.Periodic,Nothing}}
├── y: CoordinateBoundaryConditions{BoundaryCondition{Oceananigans.BoundaryConditions.Periodic,Nothing},BoundaryCondition{Oceananigans.BoundaryConditions.Periodic,Nothing}}
└── z: CoordinateBoundaryConditions{BoundaryCondition{Flux,Nothing},BoundaryCondition{Flux,Float64}}

McWilliams et al. (1997) impose a linear buoyancy gradient at the bottom along with a weak, destabilizing flux of buoyancy at the surface to faciliate spin-up from rest.

Qᵇ = 2.307e-9 # m³ s⁻², surface buoyancy flux
N² = 1.936e-5 # s⁻², initial and bottom buoyancy gradient

b_boundary_conditions = TracerBoundaryConditions(grid, top = BoundaryCondition(Flux, Qᵇ),
                                                       bottom = BoundaryCondition(Gradient, N²))
Oceananigans.FieldBoundaryConditions (NamedTuple{(:x, :y, :z)}), with boundary conditions
├── x: CoordinateBoundaryConditions{BoundaryCondition{Oceananigans.BoundaryConditions.Periodic,Nothing},BoundaryCondition{Oceananigans.BoundaryConditions.Periodic,Nothing}}
├── y: CoordinateBoundaryConditions{BoundaryCondition{Oceananigans.BoundaryConditions.Periodic,Nothing},BoundaryCondition{Oceananigans.BoundaryConditions.Periodic,Nothing}}
└── z: CoordinateBoundaryConditions{BoundaryCondition{Gradient,Float64},BoundaryCondition{Flux,Float64}}
The flux convention in Oceananigans

Note that Oceananigans uses "positive upward" conventions for all fluxes. In consequence, a negative flux at the surface drives positive velocities, and a positive flux of buoyancy drives cooling.

Coriolis parameter

McWilliams et al. (1997) use

coriolis = FPlane(f=1e-4) # s⁻¹
FPlane{Float64}: f = 1.00e-04

which is typical for mid-latitudes on Earth.

Model instantiation

We are ready to build the model. We use a fifth-order Weighted Essentially Non-Oscillatory (WENO) advection scheme and the AnisotropicMinimumDissipation model for large eddy simulation. Because our Stokes drift does not vary in $x, y$, we use UniformStokesDrift, which expects Stokes drift functions of $z, t$ only.

using Oceananigans.Advection
using Oceananigans.Buoyancy: BuoyancyTracer
using Oceananigans.StokesDrift: UniformStokesDrift

model = IncompressibleModel(
           architecture = CPU(),
              advection = WENO5(),
            timestepper = :RungeKutta3,
                   grid = grid,
                tracers = :b,
               buoyancy = BuoyancyTracer(),
               coriolis = coriolis,
                closure = AnisotropicMinimumDissipation(),
           stokes_drift = UniformStokesDrift(∂z_uˢ=∂z_uˢ),
    boundary_conditions = (u=u_boundary_conditions, b=b_boundary_conditions),
)
IncompressibleModel{CPU, Float64}(time = 0 seconds, iteration = 0) 
├── grid: RegularCartesianGrid{Float64, Periodic, Periodic, Bounded}(Nx=32, Ny=32, Nz=48)
├── tracers: (:b,)
├── closure: VerstappenAnisotropicMinimumDissipation{Float64,NamedTuple{(:b,),Tuple{Float64}},Float64,NamedTuple{(:b,),Tuple{Float64}}}
├── buoyancy: BuoyancyTracer
└── coriolis: FPlane{Float64}

Initial conditions

We make use of random noise concentrated in the upper 4 meters for buoyancy and velocity initial conditions,

Ξ(z) = randn() * exp(z / 4)

Our initial condition for buoyancy consists of a surface mixed layer 33 m deep, a deep linear stratification, plus noise,

initial_mixed_layer_depth = 33 # m
stratification(z) = z < - initial_mixed_layer_depth ? N² * z : N² * (-initial_mixed_layer_depth)

bᵢ(x, y, z) = stratification(z) + 1e-1 * Ξ(z) * N² * model.grid.Lz
bᵢ (generic function with 1 method)

The velocity initial condition in McWilliams et al. (1997) is zero Eulerian-mean velocity. This means that we must add the Stokes drift profile to the Lagrangian-mean $u$ velocity field modeled by Oceananigans.jl. We also add noise scaled by the friction velocity to $u$ and $w$.

uᵢ(x, y, z) = uˢ(z) + sqrt(abs(Qᵘ)) * 1e-1 * Ξ(z)

wᵢ(x, y, z) = sqrt(abs(Qᵘ)) * 1e-1 * Ξ(z)

set!(model, u=uᵢ, w=wᵢ, b=bᵢ)

Setting up the simulation

We use the TimeStepWizard for adaptive time-stepping with a Courant-Freidrichs-Lewy (CFL) number of 1.0,

using Oceananigans.Utils

wizard = TimeStepWizard(cfl=1.0, Δt=45.0, max_change=1.1, max_Δt=1minute)
TimeStepWizard{Float64}(1.0, Inf, 1.1, 0.5, 60.0, 0.0, 45.0)

Nice progress messaging

We define a function that prints a helpful message with maximum absolute value of $u, v, w$ and the current wall clock time.

using Oceananigans.Diagnostics, Printf

umax = FieldMaximum(abs, model.velocities.u)
vmax = FieldMaximum(abs, model.velocities.v)
wmax = FieldMaximum(abs, model.velocities.w)

wall_clock = time_ns()

function print_progress(simulation)
    model = simulation.model

    # Print a progress message
    msg = @sprintf("i: %04d, t: %s, Δt: %s, umax = (%.1e, %.1e, %.1e) ms⁻¹, wall time: %s\n",
                   model.clock.iteration,
                   prettytime(model.clock.time),
                   prettytime(wizard.Δt),
                   umax(), vmax(), wmax(),
                   prettytime(1e-9 * (time_ns() - wall_clock))
                  )

    @info msg

    return nothing
end
print_progress (generic function with 1 method)

Now we create the simulation,

simulation = Simulation(model, iteration_interval = 10,
                                               Δt = wizard,
                                        stop_time = 4hours,
                                         progress = print_progress)
Simulation{IncompressibleModel{CPU, Float64}}
├── Model clock: time = 0 seconds, iteration = 0 
├── Next time step (TimeStepWizard{Float64}): 45 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: 4 hours, stop iteration: Inf
├── Diagnostics: OrderedCollections.OrderedDict with 1 entry:
│   └── nan_checker => NaNChecker
└── Output writers: OrderedCollections.OrderedDict with no entries

Output

A field writer

We set up an output writer for the simulation that saves all velocity fields, tracer fields, and the subgrid turbulent diffusivity.

using Oceananigans.OutputWriters

output_interval = 5minutes

fields_to_output = merge(model.velocities, model.tracers, (νₑ=model.diffusivities.νₑ,))

simulation.output_writers[:fields] =
    JLD2OutputWriter(model, fields_to_output,
                     schedule = TimeInterval(output_interval),
                     prefix = "langmuir_turbulence_fields",
                     force = true)
JLD2OutputWriter scheduled on TimeInterval(5 minutes):
├── filepath: ./langmuir_turbulence_fields.jld2
├── 5 outputs: (:u, :v, :w, :b, :νₑ)
├── field slicer: FieldSlicer(:, :, :, with_halos=false)
├── array type: Array{Float32}
├── including: [:grid, :coriolis, :buoyancy, :closure]
└── max filesize: Inf YiB

An "averages" writer

We also set up output of time- and horizontally-averaged velocity field and momentum fluxes,

using Oceananigans.Fields

u, v, w = model.velocities

U = AveragedField(u, dims=(1, 2))
V = AveragedField(v, dims=(1, 2))
B = AveragedField(model.tracers.b, dims=(1, 2))

wu = AveragedField(w * u, dims=(1, 2))
wv = AveragedField(w * v, dims=(1, 2))

simulation.output_writers[:averages] =
    JLD2OutputWriter(model, (u=U, v=V, b=B, wu=wu, wv=wv),
                     schedule = AveragedTimeInterval(output_interval, window=2minutes),
                     prefix = "langmuir_turbulence_averages",
                     force = true)
JLD2OutputWriter scheduled on TimeInterval(5 minutes):
├── filepath: ./langmuir_turbulence_averages.jld2
├── 5 outputs: (:u, :v, :b, :wu, :wv) averaged on AveragedTimeInterval(window=2 minutes, stride=1, interval=5 minutes)
├── field slicer: FieldSlicer(:, :, :, with_halos=false)
├── array type: Array{Float32}
├── including: [:grid, :coriolis, :buoyancy, :closure]
└── max filesize: Inf YiB

Running the simulation

This part is easy,

run!(simulation)
[ Info: i: 0010, t: 7.250 minutes, Δt: 45 seconds, umax = (7.1e-02, 1.5e-02, 1.8e-02) ms⁻¹, wall time: 29.481 seconds
[ Info: i: 0020, t: 14.950 minutes, Δt: 49.500 seconds, umax = (7.6e-02, 2.0e-02, 2.2e-02) ms⁻¹, wall time: 31.670 seconds
[ Info: i: 0030, t: 22.625 minutes, Δt: 52.505 seconds, umax = (8.1e-02, 2.3e-02, 2.4e-02) ms⁻¹, wall time: 34.032 seconds
[ Info: i: 0040, t: 30 minutes, Δt: 49.246 seconds, umax = (8.6e-02, 2.7e-02, 2.3e-02) ms⁻¹, wall time: 36.302 seconds
[ Info: i: 0050, t: 37.339 minutes, Δt: 46.778 seconds, umax = (8.8e-02, 3.5e-02, 2.2e-02) ms⁻¹, wall time: 38.479 seconds
[ Info: i: 0060, t: 44.542 minutes, Δt: 45.418 seconds, umax = (9.2e-02, 4.0e-02, 2.5e-02) ms⁻¹, wall time: 40.671 seconds
[ Info: i: 0070, t: 51.449 minutes, Δt: 43.485 seconds, umax = (8.8e-02, 5.0e-02, 2.5e-02) ms⁻¹, wall time: 42.951 seconds
[ Info: i: 0080, t: 58.022 minutes, Δt: 45.324 seconds, umax = (8.8e-02, 5.4e-02, 3.3e-02) ms⁻¹, wall time: 45.140 seconds
[ Info: i: 0090, t: 1.083 hours, Δt: 45.564 seconds, umax = (8.6e-02, 5.1e-02, 3.3e-02) ms⁻¹, wall time: 47.426 seconds
[ Info: i: 0100, t: 1.206 hours, Δt: 46.762 seconds, umax = (8.7e-02, 5.2e-02, 3.8e-02) ms⁻¹, wall time: 49.619 seconds
[ Info: i: 0110, t: 1.327 hours, Δt: 45.985 seconds, umax = (8.7e-02, 5.7e-02, 3.4e-02) ms⁻¹, wall time: 51.819 seconds
[ Info: i: 0120, t: 1.442 hours, Δt: 46.219 seconds, umax = (8.6e-02, 5.2e-02, 3.1e-02) ms⁻¹, wall time: 54.204 seconds
[ Info: i: 0130, t: 1.564 hours, Δt: 46.251 seconds, umax = (8.7e-02, 5.6e-02, 3.6e-02) ms⁻¹, wall time: 56.402 seconds
[ Info: i: 0140, t: 1.679 hours, Δt: 46.069 seconds, umax = (8.9e-02, 5.6e-02, 3.9e-02) ms⁻¹, wall time: 58.674 seconds
[ Info: i: 0150, t: 1.800 hours, Δt: 44.843 seconds, umax = (8.6e-02, 5.6e-02, 3.6e-02) ms⁻¹, wall time: 1.014 minutes
[ Info: i: 0160, t: 1.917 hours, Δt: 46.350 seconds, umax = (9.1e-02, 5.5e-02, 3.1e-02) ms⁻¹, wall time: 1.052 minutes
[ Info: i: 0170, t: 2.037 hours, Δt: 44.020 seconds, umax = (9.2e-02, 6.2e-02, 2.8e-02) ms⁻¹, wall time: 1.088 minutes
[ Info: i: 0180, t: 2.156 hours, Δt: 43.607 seconds, umax = (9.0e-02, 6.6e-02, 2.9e-02) ms⁻¹, wall time: 1.125 minutes
[ Info: i: 0190, t: 2.275 hours, Δt: 44.422 seconds, umax = (9.1e-02, 6.7e-02, 3.0e-02) ms⁻¹, wall time: 1.162 minutes
[ Info: i: 0200, t: 2.395 hours, Δt: 44.194 seconds, umax = (8.4e-02, 7.4e-02, 3.0e-02) ms⁻¹, wall time: 1.199 minutes
[ Info: i: 0210, t: 2.513 hours, Δt: 47.379 seconds, umax = (8.7e-02, 7.2e-02, 3.7e-02) ms⁻¹, wall time: 1.236 minutes
[ Info: i: 0220, t: 2.635 hours, Δt: 46.076 seconds, umax = (8.5e-02, 7.0e-02, 2.9e-02) ms⁻¹, wall time: 1.273 minutes
[ Info: i: 0230, t: 2.750 hours, Δt: 46.950 seconds, umax = (8.8e-02, 7.1e-02, 3.2e-02) ms⁻¹, wall time: 1.310 minutes
[ Info: i: 0240, t: 2.871 hours, Δt: 45.316 seconds, umax = (8.3e-02, 7.1e-02, 2.8e-02) ms⁻¹, wall time: 1.347 minutes
[ Info: i: 0250, t: 2.997 hours, Δt: 48.140 seconds, umax = (8.4e-02, 6.8e-02, 3.0e-02) ms⁻¹, wall time: 1.383 minutes
[ Info: i: 0260, t: 3.110 hours, Δt: 47.406 seconds, umax = (8.7e-02, 6.4e-02, 3.0e-02) ms⁻¹, wall time: 1.421 minutes
[ Info: i: 0270, t: 3.230 hours, Δt: 45.888 seconds, umax = (8.9e-02, 6.7e-02, 2.8e-02) ms⁻¹, wall time: 1.457 minutes
[ Info: i: 0280, t: 3.346 hours, Δt: 44.743 seconds, umax = (8.6e-02, 6.9e-02, 3.2e-02) ms⁻¹, wall time: 1.495 minutes
[ Info: i: 0290, t: 3.469 hours, Δt: 46.702 seconds, umax = (8.4e-02, 6.8e-02, 3.2e-02) ms⁻¹, wall time: 1.533 minutes
[ Info: i: 0300, t: 3.583 hours, Δt: 47.871 seconds, umax = (8.5e-02, 6.7e-02, 2.7e-02) ms⁻¹, wall time: 1.571 minutes
[ Info: i: 0310, t: 3.706 hours, Δt: 47.165 seconds, umax = (9.0e-02, 7.3e-02, 3.1e-02) ms⁻¹, wall time: 1.607 minutes
[ Info: i: 0320, t: 3.824 hours, Δt: 44.496 seconds, umax = (9.1e-02, 7.6e-02, 2.9e-02) ms⁻¹, wall time: 1.643 minutes
[ Info: i: 0330, t: 3.941 hours, Δt: 43.921 seconds, umax = (8.7e-02, 8.1e-02, 3.6e-02) ms⁻¹, wall time: 1.681 minutes
[ Info: i: 0335, t: 4 hours, Δt: 46.120 seconds, umax = (8.4e-02, 7.9e-02, 3.7e-02) ms⁻¹, wall time: 1.700 minutes
[ Info: Simulation is stopping. Model time 4 hours has hit or exceeded simulation stop time 4 hours.

Making a neat movie

We look at the results by plotting vertical slices of $u$ and $w$, and a horizontal slice of $w$ to look for Langmuir cells.

k = searchsortedfirst(grid.zF[:], -8)

Making the coordinate arrays takes a few lines of code,

using Oceananigans.Grids

xw, yw, zw = nodes(model.velocities.w)
xu, yu, zu = nodes(model.velocities.u)

Next, we open the JLD2 file, and extract the iterations we ended up saving at,

using JLD2, Plots

fields_file = jldopen(simulation.output_writers[:fields].filepath)
averages_file = jldopen(simulation.output_writers[:averages].filepath)

iterations = parse.(Int, keys(fields_file["timeseries/t"]))
49-element Array{Int64,1}:
   0
   7
  14
  21
  27
  33
  40
  47
  54
  61
   ⋮
 279
 286
 293
 300
 307
 314
 321
 328
 335

This utility is handy for calculating nice contour intervals:

function nice_divergent_levels(c, clim; nlevels=20)
    levels = range(-clim, stop=clim, length=nlevels)
    cmax = maximum(abs, c)
    clim < cmax && (levels = vcat([-cmax], levels, [cmax]))
    return (-clim, clim), levels
end

Finally, we're ready to animate.

@info "Making an animation from the saved data..."

anim = @animate for (i, iter) in enumerate(iterations)

    @info "Drawing frame $i from iteration $iter \n"

    # Load 3D fields from fields_file
    t = fields_file["timeseries/t/$iter"]
    w_snapshot = fields_file["timeseries/w/$iter"]
    u_snapshot = fields_file["timeseries/u/$iter"]

    B_snapshot = averages_file["timeseries/b/$iter"][1, 1, :]
    U_snapshot = averages_file["timeseries/u/$iter"][1, 1, :]
    V_snapshot = averages_file["timeseries/v/$iter"][1, 1, :]
    wu_snapshot = averages_file["timeseries/wu/$iter"][1, 1, :]
    wv_snapshot = averages_file["timeseries/wv/$iter"][1, 1, :]

    # Extract slices
    wxy = w_snapshot[:, :, k]
    wxz = w_snapshot[:, 1, :]
    uxz = u_snapshot[:, 1, :]

    wlims, wlevels = nice_divergent_levels(w, 0.03)
    ulims, ulevels = nice_divergent_levels(w, 0.05)

    B_plot = plot(B_snapshot, zu,
                  label = nothing,
                  legend = :bottom,
                  xlabel = "Buoyancy (m s⁻²)",
                  ylabel = "z (m)")

    U_plot = plot([U_snapshot V_snapshot], zu,
                  label = ["\$ \\bar u \$" "\$ \\bar v \$"],
                  legend = :bottom,
                  xlabel = "Velocities (m s⁻¹)",
                  ylabel = "z (m)")

    wu_label = "\$ \\overline{wu} \$"
    wv_label = "\$ \\overline{wv} \$"

    fluxes_plot = plot([wu_snapshot, wv_snapshot], zw,
                       label = [wu_label wv_label],
                       legend = :bottom,
                       xlabel = "Momentum fluxes (m² s⁻²)",
                       ylabel = "z (m)")

    wxy_plot = contourf(xw, yw, wxy';
                        color = :balance,
                        linewidth = 0,
                        aspectratio = :equal,
                        clims = wlims,
                        levels = wlevels,
                        xlims = (0, grid.Lx),
                        ylims = (0, grid.Ly),
                        xlabel = "x (m)",
                        ylabel = "y (m)")

    wxz_plot = contourf(xw, zw, wxz';
                              color = :balance,
                          linewidth = 0,
                        aspectratio = :equal,
                              clims = wlims,
                             levels = wlevels,
                              xlims = (0, grid.Lx),
                              ylims = (-grid.Lz, 0),
                             xlabel = "x (m)",
                             ylabel = "z (m)")

    uxz_plot = contourf(xu, zu, uxz';
                              color = :balance,
                          linewidth = 0,
                        aspectratio = :equal,
                              clims = ulims,
                             levels = ulevels,
                              xlims = (0, grid.Lx),
                              ylims = (-grid.Lz, 0),
                             xlabel = "x (m)",
                             ylabel = "z (m)")

    wxy_title = @sprintf("w(x, y, t) (m s⁻¹) at z=-8 m and t = %s ", prettytime(t))
    wxz_title = @sprintf("w(x, z, t) (m s⁻¹) at y=0 m and t = %s", prettytime(t))
    uxz_title = @sprintf("u(x, z, t) (m s⁻¹) at y=0 m and t = %s", prettytime(t))

    plot(wxy_plot, B_plot, wxz_plot, U_plot, uxz_plot, fluxes_plot,
         layout = Plots.grid(3, 2, widths=(0.7, 0.3)), size = (900.5, 1000.5),
         title = [wxy_title "" wxz_title "" uxz_title ""])

    if iter == iterations[end]
        close(fields_file)
        close(averages_file)
    end
end

This page was generated using Literate.jl.