Langmuir turbulence example

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

Wagner et al., "Near-inertial waves and turbulence driven by the growth of swell", Journal of Physical Oceanography (2021)

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, Plots"
using Oceananigans
using Oceananigans.Units: minute, minutes, hours

Model set-up

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

Domain and numerical grid specification

We use a modest resolution and the same total extent as Wagner et al. 2021,

grid = RectilinearGrid(size=(32, 32, 32), extent=(128, 128, 64))
32×32×32 RectilinearGrid{Float64, Periodic, Periodic, Bounded} on CPU with 1×1×1 halo
├── Periodic x ∈ [0.0, 128.0) regularly spaced with Δx=4.0
├── Periodic y ∈ [0.0, 128.0) regularly spaced with Δy=4.0
└── Bounded  z ∈ [-64.0, 0.0] regularly spaced with Δz=2.0

The Stokes Drift profile

The surface wave Stokes drift profile prescribed in Wagner et al. 2021, 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.BuoyancyModels: 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 NonhydrostaticModel 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 and surface wave field undergoes forced 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$, Wagner et al. 2021 impose

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

u_boundary_conditions = FieldBoundaryConditions(top = FluxBoundaryCondition(Qᵘ))
Oceananigans.FieldBoundaryConditions, with boundary conditions
├── west: DefaultBoundaryCondition
├── east: DefaultBoundaryCondition
├── south: DefaultBoundaryCondition
├── north: DefaultBoundaryCondition
├── bottom: DefaultBoundaryCondition
├── top: FluxBoundaryCondition: -3.72e-5
└── immersed: FluxBoundaryCondition: Nothing

Wagner et al. 2021 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 = FieldBoundaryConditions(top = FluxBoundaryCondition(Qᵇ),
                                                bottom = GradientBoundaryCondition(N²))
Oceananigans.FieldBoundaryConditions, with boundary conditions
├── west: DefaultBoundaryCondition
├── east: DefaultBoundaryCondition
├── south: DefaultBoundaryCondition
├── north: DefaultBoundaryCondition
├── bottom: GradientBoundaryCondition: 1.936e-5
├── top: FluxBoundaryCondition: 2.307e-9
└── immersed: FluxBoundaryCondition: Nothing
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

Wagner et al. (2021) use

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

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.

model = NonhydrostaticModel(; grid,
                            advection = WENO5(),
                            timestepper = :RungeKutta3,
                            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))
NonhydrostaticModel{CPU, RectilinearGrid}(time = 0 seconds, iteration = 0)
├── grid: 32×32×32 RectilinearGrid{Float64, Periodic, Periodic, Bounded} on CPU with 3×3×3 halo
├── timestepper: RungeKutta3TimeStepper
├── tracers: b
├── closure: AnisotropicMinimumDissipation{ExplicitTimeDiscretization, NamedTuple{(:b,), Tuple{Float64}}, Float64, Nothing}
├── buoyancy: BuoyancyTracer with -ĝ = ZDirection
└── coriolis: FPlane{Float64}(f=0.0001)

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 simulation we reproduce from Wagner et al. (2021) is zero Lagrangian-mean velocity. This initial condition is consistent with a wavy, quiescent ocean suddenly impacted by winds. To this quiescent state we add noise scaled by the friction velocity to $u$ and $w$.

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

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

Setting up the simulation

simulation = Simulation(model, Δt=45.0, stop_time=4hours)
Simulation of NonhydrostaticModel{CPU, RectilinearGrid}(time = 0 seconds, iteration = 0)
├── Next time step: 45 seconds
├── Elapsed wall time: 0 seconds
├── Wall time per iteration: NaN years
├── Stop time: 4 hours
├── 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

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

wizard = TimeStepWizard(cfl=1.0, max_change=1.1, max_Δt=1minute)

simulation.callbacks[:wizard] = Callback(wizard, IterationInterval(10))
Callback of TimeStepWizard(cfl=1.0, max_Δt=60.0, min_Δt=0.0) on IterationInterval(10)

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 Printf

function print_progress(simulation)
    u, v, w = simulation.model.velocities

    # Print a progress message
    msg = @sprintf("i: %04d, t: %s, Δt: %s, umax = (%.1e, %.1e, %.1e) ms⁻¹, wall time: %s\n",
                   iteration(simulation),
                   prettytime(time(simulation)),
                   prettytime(simulation.Δt),
                   maximum(abs, u), maximum(abs, v), maximum(abs, w),
                   prettytime(simulation.run_wall_time))

    @info msg

    return nothing
end

simulation.callbacks[:progress] = Callback(print_progress, IterationInterval(10))
Callback of print_progress on IterationInterval(10)

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.

output_interval = 5minutes

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

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, νₑ)
├── 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,

u, v, w = model.velocities
b = model.tracers.b

 U = Average(u, dims=(1, 2))
 V = Average(v, dims=(1, 2))
 B = Average(b, dims=(1, 2))
wu = Average(w * u, dims=(1, 2))
wv = Average(w * v, dims=(1, 2))

simulation.output_writers[:averages] =
    JLD2OutputWriter(model, (; U, V, B, wu, 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)
├── array type: Array{Float32}
├── including: [:grid, :coriolis, :buoyancy, :closure]
└── max filesize: Inf YiB

Running the simulation

This part is easy,

run!(simulation)
[ Info: Initializing simulation...
[ Info: i: 0000, t: 0 seconds, Δt: 49.500 seconds, umax = (1.3e-03, 6.2e-04, 5.8e-04) ms⁻¹, wall time: 0 seconds
[ Info:     ... simulation initialization complete (8.915 seconds)
[ Info: Executing initial time step...
[ Info:     ... initial time step complete (1.912 minutes).
[ Info: i: 0010, t: 7.475 minutes, Δt: 54.450 seconds, umax = (2.2e-02, 1.1e-02, 2.0e-02) ms⁻¹, wall time: 2.081 minutes
[ Info: i: 0020, t: 15.907 minutes, Δt: 59.895 seconds, umax = (3.0e-02, 1.1e-02, 1.9e-02) ms⁻¹, wall time: 2.102 minutes
[ Info: i: 0030, t: 24.991 minutes, Δt: 1 minute, umax = (3.5e-02, 1.0e-02, 1.8e-02) ms⁻¹, wall time: 2.123 minutes
[ Info: i: 0040, t: 33 minutes, Δt: 1 minute, umax = (4.1e-02, 1.1e-02, 2.0e-02) ms⁻¹, wall time: 2.144 minutes
[ Info: i: 0050, t: 43 minutes, Δt: 1 minute, umax = (4.8e-02, 1.4e-02, 2.1e-02) ms⁻¹, wall time: 2.165 minutes
[ Info: i: 0060, t: 53 minutes, Δt: 1 minute, umax = (5.6e-02, 2.1e-02, 1.7e-02) ms⁻¹, wall time: 2.190 minutes
[ Info: i: 0070, t: 1.050 hours, Δt: 1 minute, umax = (6.0e-02, 2.7e-02, 2.2e-02) ms⁻¹, wall time: 2.211 minutes
[ Info: i: 0080, t: 1.217 hours, Δt: 1 minute, umax = (6.3e-02, 3.1e-02, 1.9e-02) ms⁻¹, wall time: 2.232 minutes
[ Info: i: 0090, t: 1.383 hours, Δt: 1 minute, umax = (6.4e-02, 3.4e-02, 2.3e-02) ms⁻¹, wall time: 2.252 minutes
[ Info: i: 0100, t: 1.550 hours, Δt: 1 minute, umax = (6.6e-02, 3.0e-02, 1.9e-02) ms⁻¹, wall time: 2.273 minutes
[ Info: i: 0110, t: 1.717 hours, Δt: 57.804 seconds, umax = (6.9e-02, 3.8e-02, 2.1e-02) ms⁻¹, wall time: 2.295 minutes
[ Info: i: 0120, t: 1.833 hours, Δt: 55.850 seconds, umax = (7.2e-02, 3.4e-02, 2.0e-02) ms⁻¹, wall time: 2.316 minutes
[ Info: i: 0130, t: 1.979 hours, Δt: 53.315 seconds, umax = (7.5e-02, 3.9e-02, 1.9e-02) ms⁻¹, wall time: 2.337 minutes
[ Info: i: 0140, t: 2.113 hours, Δt: 53.012 seconds, umax = (7.5e-02, 4.2e-02, 2.4e-02) ms⁻¹, wall time: 2.358 minutes
[ Info: i: 0150, t: 2.250 hours, Δt: 48.933 seconds, umax = (8.2e-02, 4.3e-02, 2.2e-02) ms⁻¹, wall time: 2.381 minutes
[ Info: i: 0160, t: 2.361 hours, Δt: 50.834 seconds, umax = (7.9e-02, 4.9e-02, 2.3e-02) ms⁻¹, wall time: 2.402 minutes
[ Info: i: 0170, t: 2.500 hours, Δt: 50.314 seconds, umax = (8.0e-02, 5.0e-02, 3.0e-02) ms⁻¹, wall time: 2.422 minutes
[ Info: i: 0180, t: 2.639 hours, Δt: 48.042 seconds, umax = (8.3e-02, 5.5e-02, 3.0e-02) ms⁻¹, wall time: 2.443 minutes
[ Info: i: 0190, t: 2.750 hours, Δt: 47.415 seconds, umax = (8.4e-02, 5.0e-02, 3.2e-02) ms⁻¹, wall time: 2.463 minutes
[ Info: i: 0200, t: 2.873 hours, Δt: 46.964 seconds, umax = (8.5e-02, 4.9e-02, 3.1e-02) ms⁻¹, wall time: 2.484 minutes
[ Info: i: 0210, t: 2.995 hours, Δt: 48.163 seconds, umax = (8.3e-02, 4.7e-02, 3.7e-02) ms⁻¹, wall time: 2.504 minutes
[ Info: i: 0220, t: 3.110 hours, Δt: 49.014 seconds, umax = (8.2e-02, 5.5e-02, 3.7e-02) ms⁻¹, wall time: 2.524 minutes
[ Info: i: 0230, t: 3.235 hours, Δt: 49.693 seconds, umax = (8.0e-02, 6.0e-02, 2.9e-02) ms⁻¹, wall time: 2.550 minutes
[ Info: i: 0240, t: 3.347 hours, Δt: 52.201 seconds, umax = (7.7e-02, 6.0e-02, 3.5e-02) ms⁻¹, wall time: 2.576 minutes
[ Info: i: 0250, t: 3.489 hours, Δt: 50.011 seconds, umax = (8.0e-02, 6.3e-02, 3.5e-02) ms⁻¹, wall time: 2.598 minutes
[ Info: i: 0260, t: 3.597 hours, Δt: 50.824 seconds, umax = (7.9e-02, 6.0e-02, 3.5e-02) ms⁻¹, wall time: 2.619 minutes
[ Info: i: 0270, t: 3.737 hours, Δt: 50.872 seconds, umax = (7.9e-02, 6.1e-02, 2.6e-02) ms⁻¹, wall time: 2.639 minutes
[ Info: i: 0280, t: 3.862 hours, Δt: 48.978 seconds, umax = (8.2e-02, 5.9e-02, 2.8e-02) ms⁻¹, wall time: 2.659 minutes
[ Info: i: 0290, t: 3.985 hours, Δt: 47.770 seconds, umax = (8.4e-02, 5.9e-02, 2.9e-02) ms⁻¹, wall time: 2.679 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 loading data from file with FieldTimeSeries, and plotting vertical slices of $u$ and $w$, and a horizontal slice of $w$ to look for Langmuir cells.

using Plots

time_series = (;
     w = FieldTimeSeries("langmuir_turbulence_fields.jld2", "w"),
     u = FieldTimeSeries("langmuir_turbulence_fields.jld2", "u"),
     B = FieldTimeSeries("langmuir_turbulence_averages.jld2", "B"),
     U = FieldTimeSeries("langmuir_turbulence_averages.jld2", "U"),
     V = FieldTimeSeries("langmuir_turbulence_averages.jld2", "V"),
    wu = FieldTimeSeries("langmuir_turbulence_averages.jld2", "wu"),
    wv = FieldTimeSeries("langmuir_turbulence_averages.jld2", "wv"))

times = time_series.w.times
xw, yw, zw = nodes(time_series.w)
xu, yu, zu = nodes(time_series.u)
WARNING: using Plots.grid in module Main conflicts with an existing identifier.

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, t) in enumerate(times)

    @info "Drawing frame $i of $(length(times))..."

    # Get fields at ith save point
     wᵢ = time_series.w[i]
     uᵢ = time_series.u[i]
     Bᵢ = time_series.B[i][1, 1, :]
     Uᵢ = time_series.U[i][1, 1, :]
     Vᵢ = time_series.V[i][1, 1, :]
    wuᵢ = time_series.wu[i][1, 1, :]
    wvᵢ = time_series.wv[i][1, 1, :]

    # Extract slices
    k = searchsortedfirst(grid.zᵃᵃᶠ[:], -8)
    wxy = wᵢ[:, :, k]
    wxz = wᵢ[:, 1, :]
    uxz = uᵢ[:, 1, :]

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

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

    U_plot = plot([Uᵢ Vᵢ], 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ᵢ, wvᵢ], 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 ""])
end

This page was generated using Literate.jl.