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 = UpwindBiasedFifthOrder(),
                            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
├── 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, y, z) = 2e-2 * randn()
shear_flow(x, y, z) = tanh(z) + noise(x, y, z)

Ri₀ = 0.1; h = 1/4
stratification(x, y, 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: ("128×1×128 Field{Face, Center, Center} on RectilinearGrid on CPU", "128×1×128 Field{Center, Face, Center} on RectilinearGrid on CPU", "128×1×129 Field{Center, Center, Face} on RectilinearGrid on CPU", "128×1×128 Field{Center, Center, Center} on RectilinearGrid on CPU", "Tuple{Int64, Int64, Int64}")

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: ("128×1×128 Field{Face, Center, Center} on RectilinearGrid on CPU", "128×1×128 Field{Center, Face, Center} on RectilinearGrid on CPU", "128×1×129 Field{Center, Center, Face} on RectilinearGrid on CPU")

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(TracerVarianceDissipationRate(model, :b))
∫χ = Integral(TracerVarianceDiffusiveTerm(model, :b))
sum! over dims (1, 2, 3) of BinaryOperation at (Center, Center, Center)
└── 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, ∫χ, ∫χᴰ)
filename = "kelvin_helmholtz"
simulation.output_writers[:nc] = NetCDFOutputWriter(model, output_fields,
                                                    filename = joinpath(@__DIR__, filename),
                                                    schedule = TimeInterval(1),
                                                    overwrite_existing = true)
NetCDFOutputWriter scheduled on TimeInterval(1 second):
├── filepath: /home/runner/work/Oceanostics.jl/Oceanostics.jl/docs/build/generated/kelvin_helmholtz.nc
├── dimensions: zC(128), zF(129), xC(128), yF(1), xF(128), yC(1), time(0)
├── 5 outputs: (Q, ∫χᴰ, Ri, b, ∫χ)
└── array type: Array{Float64}

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.276 ms,  walltime = 26.547 seconds,  walltime / timestep = 0 seconds
          |u⃗|ₘₐₓ = [1.04e+00,  0.00e+00,  3.63e-02] m/s,  advective CFL = 0.8,  diffusive CFL = 0.00019,  νₘₐₓ = 2e-05 m²/s
[ Info:     ... simulation initialization complete (20.134 seconds)
[ Info: Executing initial time step...
[ Info:     ... initial time step complete (7.786 seconds).
┌ Info: iter =    200,  [011.79%] time = 11.793 seconds,  Δt = 61.014 ms,  walltime = 1.132 minutes,  walltime / timestep = 206.749 ms
          |u⃗|ₘₐₓ = [1.02e+00,  0.00e+00,  2.67e-02] m/s,  advective CFL = 0.8,  diffusive CFL = 0.0002,  νₘₐₓ = 2e-05 m²/s
┌ Info: iter =    400,  [023.55%] time = 23.550 seconds,  Δt = 61.086 ms,  walltime = 1.406 minutes,  walltime / timestep = 82.396 ms
          |u⃗|ₘₐₓ = [1.02e+00,  0.00e+00,  4.06e-02] m/s,  advective CFL = 0.8,  diffusive CFL = 0.0002,  νₘₐₓ = 2e-05 m²/s
┌ Info: iter =    600,  [035.30%] time = 35.304 seconds,  Δt = 60.727 ms,  walltime = 1.676 minutes,  walltime / timestep = 81.044 ms
          |u⃗|ₘₐₓ = [1.02e+00,  0.00e+00,  6.32e-02] m/s,  advective CFL = 0.8,  diffusive CFL = 0.0002,  νₘₐₓ = 2e-05 m²/s
┌ Info: iter =    800,  [046.74%] time = 46.737 seconds,  Δt = 56.046 ms,  walltime = 1.937 minutes,  walltime / timestep = 78.126 ms
          |u⃗|ₘₐₓ = [1.07e+00,  0.00e+00,  1.76e-01] m/s,  advective CFL = 0.8,  diffusive CFL = 0.00018,  νₘₐₓ = 2e-05 m²/s
┌ Info: iter =   1000,  [057.00%] time = 56.997 seconds,  Δt = 49.638 ms,  walltime = 2.191 minutes,  walltime / timestep = 76.253 ms
          |u⃗|ₘₐₓ = [1.16e+00,  0.00e+00,  3.75e-01] m/s,  advective CFL = 0.8,  diffusive CFL = 0.00016,  νₘₐₓ = 2e-05 m²/s
┌ Info: iter =   1200,  [066.50%] time = 1.108 minutes,  Δt = 49.672 ms,  walltime = 2.443 minutes,  walltime / timestep = 75.623 ms
          |u⃗|ₘₐₓ = [1.15e+00,  0.00e+00,  3.80e-01] m/s,  advective CFL = 0.8,  diffusive CFL = 0.00016,  νₘₐₓ = 2e-05 m²/s
┌ Info: iter =   1400,  [076.31%] time = 1.272 minutes,  Δt = 52.477 ms,  walltime = 2.690 minutes,  walltime / timestep = 74.165 ms
          |u⃗|ₘₐₓ = [1.12e+00,  0.00e+00,  3.00e-01] m/s,  advective CFL = 0.8,  diffusive CFL = 0.00017,  νₘₐₓ = 2e-05 m²/s
┌ Info: iter =   1600,  [086.68%] time = 1.445 minutes,  Δt = 51.944 ms,  walltime = 2.935 minutes,  walltime / timestep = 73.333 ms
          |u⃗|ₘₐₓ = [1.12e+00,  0.00e+00,  3.10e-01] m/s,  advective CFL = 0.8,  diffusive CFL = 0.00017,  νₘₐₓ = 2e-05 m²/s
┌ Info: iter =   1800,  [096.49%] time = 1.608 minutes,  Δt = 49.355 ms,  walltime = 3.177 minutes,  walltime / timestep = 72.775 ms
          |u⃗|ₘₐₓ = [1.16e+00,  0.00e+00,  3.56e-01] m/s,  advective CFL = 0.8,  diffusive CFL = 0.00016,  νₘₐₓ = 2e-05 m²/s
[ Info: Simulation is stopping after running for 2.867 minutes.
[ 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)
RasterStack with dimensions: 
  Dim{:xC} Sampled{Float64} Float64[-4.9609375, -4.8828125, …, 4.8828125, 4.9609375] ForwardOrdered Regular Points,
  Dim{:yC} Sampled{Float64} Float64[1.0] ForwardOrdered Regular Points,
  Dim{:zC} Sampled{Float64} Float64[-4.9609375, -4.8828125, …, 4.8828125, 4.9609375] ForwardOrdered Regular Points,
  Ti Sampled{Float64} Float64[0.0, 1.0, …, 99.0, 100.0] ForwardOrdered Regular Points,
  Dim{:zF} Sampled{Float64} Float64[-5.0, -4.921875, …, 4.921875, 5.0] ForwardOrdered Regular Points
and 5 layers:
  :Q   Float64 dims: Dim{:xC}, Dim{:yC}, Dim{:zC}, Ti (128×1×128×101)
  :∫χᴰ Float64 dims: Ti (101)
  :Ri  Float64 dims: Dim{:xC}, Dim{:yC}, Dim{:zF}, Ti (128×1×129×101)
  :b   Float64 dims: Dim{:xC}, Dim{:yC}, Dim{:zC}, Ti (128×1×128×101)
  :∫χ  Float64 dims: Ti (101)

with 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 Julia Version 1.9.3\…
  "output time interval" => "Output was saved every 1 second."
  "date"                 => "This file was generated on 2023-10-06T16:03:38.312…
  "schedule"             => "TimeInterval"

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...);

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

n = Observable(1)

Riₙ = @lift set(ds.Ri[Ti=$n, yC=Near(0)], :xC => X, :zF => 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, yC=Near(0)], :xC => X, :zC => 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, yC=Near(0)], :xC => X, :zC => Z)
hm3 = heatmap!(ax3, bₙ; colormap = :balance, colorrange = (-2.5e-2, +2.5e-2))
Colorbar(fig[3, 3], hm3, vertical=false, height=8);

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), ds.∫χ,  label = "∫χdV")
lines!(axb, Array(times), 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 TracerVarianceDiffusiveTerm 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.