Kelvin-Helmholtz instability

This example simulates a simple 2D Kelvin-Helmholtz instability and is based on the similar Oceananigans example.

Before starting, make sure you have the required packages installed for this example, which can be done with

using Pkg
pkg"add Oceananigans, Oceanostics, CairoMakie, Rasters"

Model and simulation setup

We begin by creating a model with an isotropic diffusivity and fifth-order advection on a xz 128² grid using a buoyancy b as the active scalar. We'll work here with nondimensional quantities.

using Oceananigans

N = 128
L = 10
grid = RectilinearGrid(size=(N, N), x=(-L/2, +L/2), z=(-L/2, +L/2), topology=(Periodic, Flat, Bounded))

model = NonhydrostaticModel(grid; timestepper = :RungeKutta3,
                            advection = UpwindBiased(order=5),
                            closure = ScalarDiffusivity(ν=2e-5, κ=2e-5),
                            buoyancy = BuoyancyTracer(), tracers = :b)
NonhydrostaticModel{CPU, RectilinearGrid}(time = 0 seconds, iteration = 0)
├── grid: 128×1×128 RectilinearGrid{Float64, Periodic, Flat, Bounded} on CPU with 3×0×3 halo
├── timestepper: RungeKutta3TimeStepper
├── advection scheme: UpwindBiased(order=5)
├── tracers: b
├── closure: ScalarDiffusivity{ExplicitTimeDiscretization}(ν=2.0e-5, κ=(b=2.0e-5,))
├── buoyancy: BuoyancyTracer with ĝ = NegativeZDirection()
└── coriolis: Nothing

We use hyperbolic tangent functions for the initial conditions and set the maximum Richardson number below the threshold of 1/4. We also add some grid-scale small-amplitude noise to u to kick the instability off:

noise(x, z) = 2e-2 * randn()
shear_flow(x, z) = tanh(z) + noise(x, z)

Ri₀ = 0.1; h = 1/4
stratification(x, z) = h * Ri₀ * tanh(z / h)

set!(model, u=shear_flow, b=stratification)

Next create an adaptive-time-step simulation using the model above:

simulation = Simulation(model, Δt=0.1, stop_time=100)

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

Model diagnostics

We set-up a progress messenger using the TimedMessenger, which displays, among other information, the time step duration

using Oceanostics

progress = ProgressMessengers.TimedMessenger()
simulation.callbacks[:progress] = Callback(progress, IterationInterval(200))
Callback of Oceanostics.ProgressMessengers.TimedMessenger{Oceanostics.ProgressMessengers.AbstractProgressMessenger} on IterationInterval(200)

We can also define some useful diagnostics for of the flow, starting with the RichardsonNumber

Ri = RichardsonNumber(model)
KernelFunctionOperation at (Center, Center, Face)
├── grid: 128×1×128 RectilinearGrid{Float64, Periodic, Flat, Bounded} on CPU with 3×0×3 halo
├── kernel_function: richardson_number_ccf (generic function with 1 method)
└── arguments: ("Field", "Field", "Field", "Field", "Tuple")

We also set-up the QVelocityGradientTensorInvariant, which is usually used for visualizing vortices in the flow:

Q = QVelocityGradientTensorInvariant(model)
KernelFunctionOperation at (Center, Center, Center)
├── grid: 128×1×128 RectilinearGrid{Float64, Periodic, Flat, Bounded} on CPU with 3×0×3 halo
├── kernel_function: Q_velocity_gradient_tensor_invariant_ccc (generic function with 1 method)
└── arguments: ("Field", "Field", "Field")

Q is one of the velocity gradient tensor invariants and it measures the amount of vorticity versus the strain in the flow and, when it's positive, indicates a vortex. This method of vortex visualization is called the Q-criterion.

Let's also keep track of the amount of buoyancy mixing by measuring the buoyancy variance dissipation rate and diffusive term. When volume-integrated, these two quantities should be equal.

∫χᴰ = Integral(TracerVarianceEquation.DissipationRate(model, :b))
∫χ = Integral(TracerVarianceEquation.Diffusion(model, :b))
Integral of BinaryOperation at (Center, Center, Center) over dims (1, 2, 3)
└── operand: BinaryOperation at (Center, Center, Center)
    └── grid: 128×1×128 RectilinearGrid{Float64, Periodic, Flat, Bounded} on CPU with 3×0×3 halo

Now we write these quantities, along with b, to a NetCDF:

output_fields = (; Ri, Q, model.tracers.b, ∫χ, ∫χᴰ)

using NCDatasets
filename = "kelvin_helmholtz"
simulation.output_writers[:nc] = NetCDFWriter(model, output_fields,
                                              filename = joinpath(@__DIR__, filename),
                                              schedule = TimeInterval(1),
                                              overwrite_existing = true)
NetCDFWriter scheduled on TimeInterval(1 second):
├── filepath: build/generated/kelvin_helmholtz.nc
├── dimensions: time(0), x_faa(128), x_caa(128), z_aaf(129), z_aac(128)
├── 5 outputs: (Q, ∫χᴰ, Ri, b, ∫χ)
├── array_type: Array{Float32}
├── file_splitting: NoFileSplitting
└── file size: 29.2 KiB

Run the simulation and process results

To run the simulation:

run!(simulation)
[ Info: Initializing simulation...
┌ Info: iter =      0,  [000.00%] time = 0 seconds,  Δt = 59.067 ms,  walltime = 2.467 minutes,  walltime / timestep = 0 seconds
          |u⃗|ₘₐₓ = [1.04e+00,  0.00e+00,  3.08e-02] m/s,  advective CFL = 0.8,  diffusive CFL = 0.00019,  νₘₐₓ = 2e-05 m²/s
[ Info:     ... simulation initialization complete (18.059 seconds)
[ Info: Executing initial time step...
[ Info:     ... initial time step complete (3.130 seconds).
┌ Info: iter =    200,  [011.80%] time = 11.796 seconds,  Δt = 61.218 ms,  walltime = 2.828 minutes,  walltime / timestep = 108.351 ms
          |u⃗|ₘₐₓ = [1.02e+00,  0.00e+00,  2.86e-02] m/s,  advective CFL = 0.8,  diffusive CFL = 0.0002,  νₘₐₓ = 2e-05 m²/s
┌ Info: iter =    400,  [023.55%] time = 23.552 seconds,  Δt = 61.301 ms,  walltime = 2.867 minutes,  walltime / timestep = 11.645 ms
          |u⃗|ₘₐₓ = [1.02e+00,  0.00e+00,  2.64e-02] m/s,  advective CFL = 0.8,  diffusive CFL = 0.0002,  νₘₐₓ = 2e-05 m²/s
┌ Info: iter =    600,  [035.31%] time = 35.307 seconds,  Δt = 61.401 ms,  walltime = 2.904 minutes,  walltime / timestep = 11.110 ms
          |u⃗|ₘₐₓ = [1.02e+00,  0.00e+00,  2.58e-02] m/s,  advective CFL = 0.8,  diffusive CFL = 0.0002,  νₘₐₓ = 2e-05 m²/s
┌ Info: iter =    800,  [047.06%] time = 47.061 seconds,  Δt = 61.027 ms,  walltime = 2.939 minutes,  walltime / timestep = 10.457 ms
          |u⃗|ₘₐₓ = [1.02e+00,  0.00e+00,  4.38e-02] m/s,  advective CFL = 0.8,  diffusive CFL = 0.0002,  νₘₐₓ = 2e-05 m²/s
┌ Info: iter =   1000,  [058.82%] time = 58.824 seconds,  Δt = 58.728 ms,  walltime = 2.972 minutes,  walltime / timestep = 9.948 ms
          |u⃗|ₘₐₓ = [1.04e+00,  0.00e+00,  1.18e-01] m/s,  advective CFL = 0.8,  diffusive CFL = 0.00019,  νₘₐₓ = 2e-05 m²/s
┌ Info: iter =   1200,  [069.62%] time = 1.160 minutes,  Δt = 51.729 ms,  walltime = 3.004 minutes,  walltime / timestep = 9.642 ms
          |u⃗|ₘₐₓ = [1.13e+00,  0.00e+00,  3.22e-01] m/s,  advective CFL = 0.8,  diffusive CFL = 0.00017,  νₘₐₓ = 2e-05 m²/s
┌ Info: iter =   1400,  [079.24%] time = 1.321 minutes,  Δt = 48.733 ms,  walltime = 3.036 minutes,  walltime / timestep = 9.458 ms
          |u⃗|ₘₐₓ = [1.17e+00,  0.00e+00,  3.91e-01] m/s,  advective CFL = 0.8,  diffusive CFL = 0.00016,  νₘₐₓ = 2e-05 m²/s
┌ Info: iter =   1600,  [088.97%] time = 1.483 minutes,  Δt = 51.356 ms,  walltime = 3.066 minutes,  walltime / timestep = 9.044 ms
          |u⃗|ₘₐₓ = [1.13e+00,  0.00e+00,  3.21e-01] m/s,  advective CFL = 0.8,  diffusive CFL = 0.00017,  νₘₐₓ = 2e-05 m²/s
┌ Info: iter =   1800,  [099.26%] time = 1.654 minutes,  Δt = 52.659 ms,  walltime = 3.098 minutes,  walltime / timestep = 9.762 ms
          |u⃗|ₘₐₓ = [1.11e+00,  0.00e+00,  2.91e-01] m/s,  advective CFL = 0.8,  diffusive CFL = 0.00017,  νₘₐₓ = 2e-05 m²/s
[ Info: Simulation is stopping after running for 40.214 seconds.
[ Info: Simulation time 1.667 minutes equals or exceeds stop time 1.667 minutes.

Now we'll read the results using Rasters.jl, which works somewhat similarly to Python's Xarray and can speed-up the work the workflow

using Rasters

ds = RasterStack(simulation.output_writers[:nc].filepath)
128×128×129×101×128 RasterStack
├─────────────────────────────────┴────────────────────────────────────── dims ┐
  ↓ x_caa Sampled{Float64} [-4.9609375, …, 4.9609375] ForwardOrdered Regular Points,
  → z_aac Sampled{Float64} [-4.9609375, …, 4.9609375] ForwardOrdered Regular Points,
  ↗ z_aaf Sampled{Float64} [-5.0, …, 5.0] ForwardOrdered Regular Points,
  ⬔ Ti Sampled{Float64} [0.0, …, 100.0] ForwardOrdered Regular Points,
  ◩ x_faa Sampled{Float64} [-5.0, …, 4.921875] ForwardOrdered Regular Points
├────────────────────────────────────────────────────────────────────── layers ┤
  :Δx_caa eltype: Union{Missing, Float32} dims: x_caa size: 128
  :Δx_faa eltype: Union{Missing, Float32} dims: x_faa size: 128
  :Δz_aac eltype: Union{Missing, Float32} dims: z_aac size: 128
  :Δz_aaf eltype: Union{Missing, Float32} dims: z_aaf size: 129
  :Q      eltype: Union{Missing, Float32} dims: x_caa, z_aac, Ti size: 128×128×101
  :Ri     eltype: Union{Missing, Float32} dims: x_caa, z_aaf, Ti size: 128×129×101
  :b      eltype: Union{Missing, Float32} dims: x_caa, z_aac, Ti size: 128×128×101
  :∫χ     eltype: Union{Missing, Float32} dims: Ti size: 101
  :∫χᴰ    eltype: Union{Missing, Float32} dims: Ti size: 101
├──────────────────────────────────────────────────────────────────── metadata ┤
  Metadata{Rasters.NCDsource} of Dict{String, Any} with 6 entries:
  "interval"             => 1.0
  "Oceananigans"         => "This file was generated using "
  "Julia"                => "This file was generated using "
  "output time interval" => "Output was saved every 1 second."
  "date"                 => "This file was generated on 2026-01-19T13:54:11.523…
  "schedule"             => "TimeInterval"
├────────────────────────────────────────────────────────────────────── raster ┤
  missingval: missing
  extent: Extent(x_caa = (-4.9609375, 4.9609375), z_aac = (-4.9609375, 4.9609375), z_aaf = (-5.0, 5.0), Ti = (0.0, 100.0), x_faa = (-5.0, 4.921875))
└──────────────────────────────────────────────────────────────────────────────┘

We now use Makie to create the figure and its axes

using CairoMakie

set_theme!(Theme(fontsize = 24))
fig = Figure()

kwargs = (xlabel="x", ylabel="z", height=150, width=250)
ax1 = Axis(fig[2, 1]; title = "Ri", kwargs...)
ax2 = Axis(fig[2, 2]; title = "Q", kwargs...)
ax3 = Axis(fig[2, 3]; title = "b", kwargs...);
Precompiling packages...
    511.3 msFilePathsGlobExt (serial)
  1 dependency successfully precompiled in 1 seconds
Precompiling packages...
    417.3 msDistancesChainRulesCoreExt (serial)
  1 dependency successfully precompiled in 0 seconds
Precompiling packages...
    984.2 msTaylorSeriesIAExt (serial)
  1 dependency successfully precompiled in 1 seconds
Precompiling packages...
   5797.1 msOceananigansMakieExt (serial)
  1 dependency successfully precompiled in 6 seconds

Next we use Observables to lift the values and plot heatmaps and their colorbars

n = Observable(1)

Riₙ = @lift set(ds.Ri[Ti=$n, y_aca=Near(0)], :x_caa => X, :z_aaf => Z)
hm1 = heatmap!(ax1, Riₙ; colormap = :bwr, colorrange = (-1, +1))
Colorbar(fig[3, 1], hm1, vertical=false, height=8)

Qₙ = @lift set(ds.Q[Ti=$n, y_aca=Near(0)], :x_caa => X, :z_aac => Z)
hm2 = heatmap!(ax2, Qₙ; colormap = :inferno, colorrange = (0, 0.2))
Colorbar(fig[3, 2], hm2, vertical=false, height=8)

bₙ = @lift set(ds.b[Ti=$n, y_aca=Near(0)], :x_caa => X, :z_aac => Z)
hm3 = heatmap!(ax3, bₙ; colormap = :balance, colorrange = (-2.5e-2, +2.5e-2))
Colorbar(fig[3, 3], hm3, vertical=false, height=8);
┌ Warning: (DimensionalData.Dimensions.Dim{:y_aca},) dims were not found in object.
@ DimensionalData.Dimensions ~/.julia/packages/DimensionalData/FWnw9/src/Dimensions/primitives.jl:852
┌ Warning: (DimensionalData.Dimensions.Dim{:y_aca},) dims were not found in object.
@ DimensionalData.Dimensions ~/.julia/packages/DimensionalData/FWnw9/src/Dimensions/primitives.jl:852
┌ Warning: (DimensionalData.Dimensions.Dim{:y_aca},) dims were not found in object.
@ DimensionalData.Dimensions ~/.julia/packages/DimensionalData/FWnw9/src/Dimensions/primitives.jl:852

We now plot the time evolution of our integrated quantities

axb = Axis(fig[4, 1:3]; xlabel="Time", height=100)
times = dims(ds, :Ti)
lines!(axb, Array(times), Array(ds.∫χ),  label = "∫χdV")
lines!(axb, Array(times), Array(ds.∫χᴰ), label = "∫χᴰdV", linestyle=:dash)
axislegend(position=:lb, labelsize=14)
Makie.Legend()

Now we mark the time by placing a vertical line in the bottom panel and adding a helpful title

tₙ = @lift times[$n]
vlines!(axb, tₙ, color=:black, linestyle=:dash)

title = @lift "Time = " * string(round(times[$n], digits=2))
fig[1, 1:3] = Label(fig, title, fontsize=24, tellwidth=false);

Finally, we adjust the figure dimensions to fit all the panels and record a movie

resize_to_layout!(fig)

@info "Animating..."
record(fig, filename * ".mp4", 1:length(times), framerate=10) do i
       n[] = i
end
"kelvin_helmholtz.mp4"

Similarly to the kinetic energy dissipation rate (see the Two-dimensional turbulence example), TracerVarianceDissipationRate and TracerVarianceDiffusion are implemented with a energy-conserving formulation, which means that (for NoFlux boundary conditions) their volume-integral should be exactly (up to machine precision) the same.


This page was generated using Literate.jl.