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, CairoMakie"
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 3×3×3 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.BuoyancyFormulations: 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, include GPU() in the RectilinearGrid constructor above.

The Stokes drift profile is

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

and its z-derivative is

∂z_uˢ(z, t) = 1 / vertical_scale * Uˢ * exp(z / vertical_scale)
∂z_uˢ (generic function with 1 method)
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.

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 $z = 0$, Wagner et al. (2021) impose

τx = -3.72e-5 # m² s⁻², surface kinematic momentum flux

u_boundary_conditions = FieldBoundaryConditions(top = FluxBoundaryCondition(τx))
Oceananigans.FieldBoundaryConditions, with boundary conditions
├── west: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── east: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── south: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── north: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── bottom: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── top: FluxBoundaryCondition: -3.72e-5
└── immersed: DefaultBoundaryCondition (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.

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

b_boundary_conditions = FieldBoundaryConditions(top = FluxBoundaryCondition(Jᵇ),
                                                bottom = GradientBoundaryCondition(N²))
Oceananigans.FieldBoundaryConditions, with boundary conditions
├── west: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── east: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── south: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── north: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── bottom: GradientBoundaryCondition: 1.936e-5
├── top: FluxBoundaryCondition: 2.307e-8
└── immersed: DefaultBoundaryCondition (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, coriolis,
                            advection = WENO(),
                            timestepper = :RungeKutta3,
                            tracers = :b,
                            buoyancy = BuoyancyTracer(),
                            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
├── advection scheme: WENO(order=5)
├── tracers: b
├── closure: AnisotropicMinimumDissipation{ExplicitTimeDiscretization, @NamedTuple{b::Float64}, Float64, Nothing}
├── buoyancy: BuoyancyTracer with ĝ = NegativeZDirection()
└── 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(τx))
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 days
├── Stop time: 4 hours
├── 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

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

conjure_time_step_wizard!(simulation, cfl=1.0, max_Δt=1minute)

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 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(progress, IterationInterval(20))
Callback of progress on IterationInterval(20)

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),
                     filename = "langmuir_turbulence_fields.jld2",
                     overwrite_existing = true)
JLD2OutputWriter scheduled on TimeInterval(5 minutes):
├── filepath: langmuir_turbulence_fields.jld2
├── 5 outputs: (u, v, w, b, νₑ)
├── array type: Array{Float64}
├── including: [:grid, :coriolis, :buoyancy, :closure]
├── file_splitting: NoFileSplitting
└── file size: 48.0 KiB

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),
                     filename = "langmuir_turbulence_averages.jld2",
                     overwrite_existing = 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{Float64}
├── including: [:grid, :coriolis, :buoyancy, :closure]
├── file_splitting: NoFileSplitting
└── file size: 42.7 KiB

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, 5.2e-04, 6.1e-04) ms⁻¹, wall time: 0 seconds
[ Info:     ... simulation initialization complete (14.984 seconds)
[ Info: Executing initial time step...
[ Info:     ... initial time step complete (6.273 seconds).
[ Info: i: 0020, t: 15.907 minutes, Δt: 59.895 seconds, umax = (2.3e-02, 8.6e-03, 1.9e-02) ms⁻¹, wall time: 24.104 seconds
[ Info: i: 0040, t: 33 minutes, Δt: 1 minute, umax = (3.5e-02, 8.8e-03, 1.6e-02) ms⁻¹, wall time: 24.972 seconds
[ Info: i: 0060, t: 53 minutes, Δt: 1 minute, umax = (4.8e-02, 1.3e-02, 1.6e-02) ms⁻¹, wall time: 25.812 seconds
[ Info: i: 0080, t: 1.183 hours, Δt: 54.549 seconds, umax = (5.5e-02, 1.8e-02, 1.5e-02) ms⁻¹, wall time: 26.654 seconds
[ Info: i: 0100, t: 1.461 hours, Δt: 50.690 seconds, umax = (6.0e-02, 2.2e-02, 1.9e-02) ms⁻¹, wall time: 27.486 seconds
[ Info: i: 0120, t: 1.722 hours, Δt: 48.840 seconds, umax = (6.3e-02, 2.8e-02, 2.2e-02) ms⁻¹, wall time: 28.343 seconds
[ Info: i: 0140, t: 1.956 hours, Δt: 45.883 seconds, umax = (6.4e-02, 3.5e-02, 2.4e-02) ms⁻¹, wall time: 29.318 seconds
[ Info: i: 0160, t: 2.191 hours, Δt: 43.470 seconds, umax = (6.5e-02, 4.1e-02, 2.8e-02) ms⁻¹, wall time: 30.687 seconds
[ Info: i: 0180, t: 2.417 hours, Δt: 41.328 seconds, umax = (7.1e-02, 4.1e-02, 3.1e-02) ms⁻¹, wall time: 31.550 seconds
[ Info: i: 0200, t: 2.628 hours, Δt: 38.679 seconds, umax = (7.1e-02, 4.3e-02, 3.1e-02) ms⁻¹, wall time: 32.385 seconds
[ Info: i: 0220, t: 2.821 hours, Δt: 35.825 seconds, umax = (6.7e-02, 4.8e-02, 3.3e-02) ms⁻¹, wall time: 33.141 seconds
[ Info: i: 0240, t: 3 hours, Δt: 36.805 seconds, umax = (7.0e-02, 5.0e-02, 3.1e-02) ms⁻¹, wall time: 33.930 seconds
[ Info: i: 0260, t: 3.198 hours, Δt: 37.877 seconds, umax = (7.1e-02, 5.1e-02, 3.3e-02) ms⁻¹, wall time: 34.771 seconds
[ Info: i: 0280, t: 3.406 hours, Δt: 39.381 seconds, umax = (7.0e-02, 5.1e-02, 3.2e-02) ms⁻¹, wall time: 35.582 seconds
[ Info: i: 0300, t: 3.616 hours, Δt: 38.615 seconds, umax = (7.1e-02, 5.5e-02, 3.3e-02) ms⁻¹, wall time: 36.340 seconds
[ Info: i: 0320, t: 3.826 hours, Δt: 37.194 seconds, umax = (7.0e-02, 5.0e-02, 3.0e-02) ms⁻¹, wall time: 37.144 seconds
[ Info: Simulation is stopping after running for 37.856 seconds.
[ Info: Simulation time 4 hours equals or exceeds 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 CairoMakie

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

We are now ready to animate using Makie. We use Makie's Observable to animate the data. To dive into how Observables work we refer to Makie.jl's Documentation.

n = Observable(1)

wxy_title = @lift string("w(x, y, t) at z=-8 m and t = ", prettytime(times[$n]))
wxz_title = @lift string("w(x, z, t) at y=0 m and t = ", prettytime(times[$n]))
uxz_title = @lift string("u(x, z, t) at y=0 m and t = ", prettytime(times[$n]))

fig = Figure(size = (850, 850))

ax_B = Axis(fig[1, 4];
            xlabel = "Buoyancy (m s⁻²)",
            ylabel = "z (m)")

ax_U = Axis(fig[2, 4];
            xlabel = "Velocities (m s⁻¹)",
            ylabel = "z (m)",
            limits = ((-0.07, 0.07), nothing))

ax_fluxes = Axis(fig[3, 4];
                 xlabel = "Momentum fluxes (m² s⁻²)",
                 ylabel = "z (m)",
                 limits = ((-3.5e-5, 3.5e-5), nothing))

ax_wxy = Axis(fig[1, 1:2];
              xlabel = "x (m)",
              ylabel = "y (m)",
              aspect = DataAspect(),
              limits = ((0, grid.Lx), (0, grid.Ly)),
              title = wxy_title)

ax_wxz = Axis(fig[2, 1:2];
              xlabel = "x (m)",
              ylabel = "z (m)",
              aspect = AxisAspect(2),
              limits = ((0, grid.Lx), (-grid.Lz, 0)),
              title = wxz_title)

ax_uxz = Axis(fig[3, 1:2];
              xlabel = "x (m)",
              ylabel = "z (m)",
              aspect = AxisAspect(2),
              limits = ((0, grid.Lx), (-grid.Lz, 0)),
              title = uxz_title)


wₙ = @lift time_series.w[$n]
uₙ = @lift time_series.u[$n]
Bₙ = @lift view(time_series.B[$n], 1, 1, :)
Uₙ = @lift view(time_series.U[$n], 1, 1, :)
Vₙ = @lift view(time_series.V[$n], 1, 1, :)
wuₙ = @lift view(time_series.wu[$n], 1, 1, :)
wvₙ = @lift view(time_series.wv[$n], 1, 1, :)

k = searchsortedfirst(znodes(grid, Face(); with_halos=true), -8)
wxyₙ = @lift view(time_series.w[$n], :, :, k)
wxzₙ = @lift view(time_series.w[$n], :, 1, :)
uxzₙ = @lift view(time_series.u[$n], :, 1, :)

wlims = (-0.03, 0.03)
ulims = (-0.05, 0.05)

lines!(ax_B, Bₙ)

lines!(ax_U, Uₙ; label = L"\bar{u}")
lines!(ax_U, Vₙ; label = L"\bar{v}")
axislegend(ax_U; position = :rb)

lines!(ax_fluxes, wuₙ; label = L"mean $wu$")
lines!(ax_fluxes, wvₙ; label = L"mean $wv$")
axislegend(ax_fluxes; position = :rb)

hm_wxy = heatmap!(ax_wxy, wxyₙ;
                  colorrange = wlims,
                  colormap = :balance)

Colorbar(fig[1, 3], hm_wxy; label = "m s⁻¹")

hm_wxz = heatmap!(ax_wxz, wxzₙ;
                  colorrange = wlims,
                  colormap = :balance)

Colorbar(fig[2, 3], hm_wxz; label = "m s⁻¹")

ax_uxz = heatmap!(ax_uxz, uxzₙ;
                  colorrange = ulims,
                  colormap = :balance)

Colorbar(fig[3, 3], ax_uxz; label = "m s⁻¹")

fig

And, finally, we record a movie.

frames = 1:length(times)

record(fig, "langmuir_turbulence.mp4", frames, framerate=8) do i
    n[] = i
end


This page was generated using Literate.jl.