Horizontal convection example

In "horizontal convection", a non-uniform buoyancy is imposed on top of an initially resting fluid.

This example demonstrates:

  • How to use computed Fields for output.
  • How to post-process saved output using FieldTimeSeries.

Install dependencies

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

using Pkg
pkg"add Oceananigans, JLD2, Plots"

Horizontal convection

We consider here two-dimensional horizontal convection of an incompressible flow $\boldsymbol{u} = (u, w)$ on the $(x, z)$-plane ($-L_x/2 \le x \le L_x/2$ and $-H \le z \le 0$). The flow evolves under the effect of gravity. The only forcing on the fluid comes from a prescribed, non-uniform buoyancy at the top-surface of the domain.

We start by importing Oceananigans and Printf.

using Oceananigans
using Printf

The grid

H = 1.0          # vertical domain extent
Lx = 2H          # horizontal domain extent
Nx, Nz = 128, 64 # horizontal, vertical resolution

grid = RectilinearGrid(size = (Nx, Nz),
                          x = (-Lx/2, Lx/2),
                          z = (-H, 0),
                       halo = (3, 3),
                   topology = (Bounded, Flat, Bounded))
128×1×64 RectilinearGrid{Float64, Bounded, Flat, Bounded} on CPU with 3×0×3 halo
├── Bounded  x ∈ [-1.0, 1.0] regularly spaced with Δx=0.015625
├── Flat y
└── Bounded  z ∈ [-1.0, 0.0] regularly spaced with Δz=0.015625

Boundary conditions

At the surface, the imposed buoyancy is

\[b(x, z = 0, t) = - b_* \cos (2 \pi x / L_x) \, ,\]

while zero-flux boundary conditions are imposed on all other boundaries. We use free-slip boundary conditions on $u$ and $w$ everywhere.

b★ = 1.0

@inline bˢ(x, y, t, p) = - p.b★ * cos(2π * x / p.Lx)

b_bcs = FieldBoundaryConditions(top = ValueBoundaryCondition(bˢ, parameters=(; b★, Lx)))
Oceananigans.FieldBoundaryConditions, with boundary conditions
├── west: DefaultBoundaryCondition
├── east: DefaultBoundaryCondition
├── south: DefaultBoundaryCondition
├── north: DefaultBoundaryCondition
├── bottom: DefaultBoundaryCondition
├── top: ValueBoundaryCondition: ContinuousBoundaryFunction bˢ at (Nothing, Nothing, Nothing)
└── immersed: FluxBoundaryCondition: Nothing

Non-dimensional control parameters and Turbulence closure

The problem is characterized by three non-dimensional parameters. The first is the domain's aspect ratio, $L_x / H$ and the other two are the Rayleigh ($Ra$) and Prandtl ($Pr$) numbers:

\[Ra = \frac{b_* L_x^3}{\nu \kappa} \, , \quad \text{and}\, \quad Pr = \frac{\nu}{\kappa} \, .\]

The Prandtl number expresses the ratio of momentum over heat diffusion; the Rayleigh number is a measure of the relative importance of gravity over viscosity in the momentum equation.

For a domain with a given extent, the nondimensional values of $Ra$ and $Pr$ uniquely determine the viscosity and diffusivity, i.e.,

\[\nu = \sqrt{\frac{Pr b_* L_x^3}{Ra}} \quad \text{and} \quad \kappa = \sqrt{\frac{b_* L_x^3}{Pr Ra}} \, .\]

We use isotropic viscosity and diffusivities, ν and κ whose values are obtain from the prescribed $Ra$ and $Pr$ numbers. Here, we use $Pr = 1$ and $Ra = 10^8$:

Pr = 1.0    # Prandtl number
Ra = 1e8    # Rayleigh number

ν = sqrt(Pr * b★ * Lx^3 / Ra)  # Laplacian viscosity
κ = ν * Pr                     # Laplacian diffusivity

Model instantiation

We instantiate the model with the fifth-order WENO advection scheme, a 3rd order Runge-Kutta time-stepping scheme, and a BuoyancyTracer.

model = NonhydrostaticModel(; grid,
              advection = WENO5(),
            timestepper = :RungeKutta3,
                tracers = :b,
               buoyancy = BuoyancyTracer(),
                closure = ScalarDiffusivity(ν=ν, κ=κ),
    boundary_conditions = (; b=b_bcs))
NonhydrostaticModel{CPU, RectilinearGrid}(time = 0 seconds, iteration = 0)
├── grid: 128×1×64 RectilinearGrid{Float64, Bounded, Flat, Bounded} on CPU with 3×0×3 halo
├── timestepper: RungeKutta3TimeStepper
├── tracers: b
├── closure: ScalarDiffusivity{ExplicitTimeDiscretization}(ν=0.000282843, κ=(b=0.000282843,))
├── buoyancy: BuoyancyTracer with -ĝ = ZDirection
└── coriolis: Nothing

Simulation set-up

We set up a simulation that runs up to $t = 40$ with a JLD2OutputWriter that saves the flow speed, $\sqrt{u^2 + w^2}$, the buoyancy, $b$, andthe vorticity, $\partial_z u - \partial_x w$.

simulation = Simulation(model, Δt=1e-2, stop_time=40.0)
Simulation of NonhydrostaticModel{CPU, RectilinearGrid}(time = 0 seconds, iteration = 0)
├── Next time step: 10 ms
├── Elapsed wall time: 0 seconds
├── Wall time per iteration: NaN years
├── Stop time: 40 seconds
├── 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

The TimeStepWizard

The TimeStepWizard manages the time-step adaptively, keeping the Courant-Freidrichs-Lewy (CFL) number close to 0.75 while ensuring the time-step does not increase beyond the maximum allowable value for numerical stability.

wizard = TimeStepWizard(cfl=0.75, max_change=1.2, max_Δt=1e-1)

simulation.callbacks[:wizard] = Callback(wizard, IterationInterval(50))
Callback of TimeStepWizard(cfl=0.75, max_Δt=0.1, min_Δt=0.0) on IterationInterval(50)

A progress messenger

We write a function that prints out a helpful progress message while the simulation runs.

progress(sim) = @printf("i: % 6d, sim time: % 1.3f, wall time: % 10s, Δt: % 1.4f, CFL: %.2e\n",
                        iteration(sim), time(sim), prettytime(sim.run_wall_time),
                        sim.Δt, AdvectiveCFL(sim.Δt)(sim.model))

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

Output

We use computed Fields to diagnose and output the total flow speed, the vorticity, $\zeta$, and the buoyancy, $b$. Note that computed Fields take "AbstractOperations" on Fields as input:

u, v, w = model.velocities # unpack velocity `Field`s
b = model.tracers.b        # unpack buoyancy `Field`

# total flow speed
s = sqrt(u^2 + w^2)

# y-component of vorticity
ζ = ∂z(u) - ∂x(w)

We create a JLD2OutputWriter that saves the speed, and the vorticity. We then add the JLD2OutputWriter to the simulation.

saved_output_prefix = "horizontal_convection"
saved_output_filename = saved_output_prefix * ".jld2"

simulation.output_writers[:fields] = JLD2OutputWriter(model, (; s, b, ζ),
                                                      schedule = TimeInterval(0.5),
                                                      prefix = saved_output_prefix,
                                                      force = true)

Ready to press the big red button:

run!(simulation)
[ Info: Initializing simulation...
i:      0, sim time:  0.000, wall time:  0 seconds, Δt:  0.0120, CFL: 0.00e+00
[ Info:     ... simulation initialization complete (12.053 seconds)
[ Info: Executing initial time step...
[ Info:     ... initial time step complete (1.318 minutes).
i:     50, sim time:  0.596, wall time: 1.564 minutes, Δt:  0.0144, CFL: 6.08e-03
i:    100, sim time:  1.302, wall time: 1.608 minutes, Δt:  0.0173, CFL: 2.63e-02
i:    150, sim time:  2.156, wall time: 1.653 minutes, Δt:  0.0207, CFL: 6.94e-02
i:    200, sim time:  3.166, wall time: 1.698 minutes, Δt:  0.0249, CFL: 1.52e-01
i:    250, sim time:  4.373, wall time: 1.742 minutes, Δt:  0.0299, CFL: 3.54e-01
i:    300, sim time:  5.828, wall time: 1.786 minutes, Δt:  0.0193, CFL: 7.50e-01
i:    350, sim time:  6.789, wall time: 1.830 minutes, Δt:  0.0172, CFL: 7.50e-01
i:    400, sim time:  7.638, wall time: 1.874 minutes, Δt:  0.0171, CFL: 7.50e-01
i:    450, sim time:  8.478, wall time: 1.918 minutes, Δt:  0.0167, CFL: 7.50e-01
i:    500, sim time:  9.300, wall time: 1.964 minutes, Δt:  0.0179, CFL: 7.50e-01
i:    550, sim time:  10.179, wall time: 2.008 minutes, Δt:  0.0166, CFL: 7.50e-01
i:    600, sim time:  10.999, wall time: 2.052 minutes, Δt:  0.0171, CFL: 7.50e-01
i:    650, sim time:  11.825, wall time: 2.097 minutes, Δt:  0.0190, CFL: 7.50e-01
i:    700, sim time:  12.729, wall time: 2.141 minutes, Δt:  0.0191, CFL: 7.50e-01
i:    750, sim time:  13.653, wall time: 2.193 minutes, Δt:  0.0197, CFL: 7.50e-01
i:    800, sim time:  14.618, wall time: 2.241 minutes, Δt:  0.0205, CFL: 7.50e-01
i:    850, sim time:  15.603, wall time: 2.285 minutes, Δt:  0.0214, CFL: 7.50e-01
i:    900, sim time:  16.650, wall time: 2.334 minutes, Δt:  0.0212, CFL: 7.50e-01
i:    950, sim time:  17.691, wall time: 2.385 minutes, Δt:  0.0213, CFL: 7.50e-01
i:   1000, sim time:  18.713, wall time: 2.428 minutes, Δt:  0.0220, CFL: 7.50e-01
i:   1050, sim time:  19.786, wall time: 2.472 minutes, Δt:  0.0233, CFL: 7.50e-01
i:   1100, sim time:  20.919, wall time: 2.518 minutes, Δt:  0.0246, CFL: 7.50e-01
i:   1150, sim time:  22.098, wall time: 2.573 minutes, Δt:  0.0214, CFL: 7.50e-01
i:   1200, sim time:  23.150, wall time: 2.638 minutes, Δt:  0.0211, CFL: 7.50e-01
i:   1250, sim time:  24.190, wall time: 2.707 minutes, Δt:  0.0251, CFL: 7.50e-01
i:   1300, sim time:  25.427, wall time: 2.774 minutes, Δt:  0.0261, CFL: 7.50e-01
i:   1350, sim time:  26.682, wall time: 2.841 minutes, Δt:  0.0294, CFL: 7.50e-01
i:   1400, sim time:  28.088, wall time: 2.894 minutes, Δt:  0.0220, CFL: 7.50e-01
i:   1450, sim time:  29.176, wall time: 2.941 minutes, Δt:  0.0188, CFL: 7.50e-01
i:   1500, sim time:  30.094, wall time: 2.987 minutes, Δt:  0.0211, CFL: 7.50e-01
i:   1550, sim time:  31.106, wall time: 3.035 minutes, Δt:  0.0206, CFL: 7.50e-01
i:   1600, sim time:  32.103, wall time: 3.124 minutes, Δt:  0.0215, CFL: 7.50e-01
i:   1650, sim time:  33.129, wall time: 3.176 minutes, Δt:  0.0235, CFL: 7.50e-01
i:   1700, sim time:  34.258, wall time: 3.221 minutes, Δt:  0.0225, CFL: 7.50e-01
i:   1750, sim time:  35.361, wall time: 3.265 minutes, Δt:  0.0250, CFL: 7.50e-01
i:   1800, sim time:  36.550, wall time: 3.312 minutes, Δt:  0.0295, CFL: 7.50e-01
i:   1850, sim time:  38.000, wall time: 3.356 minutes, Δt:  0.0284, CFL: 7.50e-01
i:   1900, sim time:  39.397, wall time: 3.400 minutes, Δt:  0.0300, CFL: 7.50e-01
[ Info: Simulation is stopping. Model time 40 seconds has hit or exceeded simulation stop time 40 seconds.

Load saved output, process, visualize

We animate the results by opening the JLD2 file, extracting data for the iterations we ended up saving at, and ploting the saved fields. From the saved buoyancy field we compute the buoyancy dissipation, $\chi = \kappa |\boldsymbol{\nabla} b|^2$, and plot that also.

To start we load the saved fields are FieldTimeSeries and prepare for animating the flow by creating coordinate arrays that each field lives on.

using JLD2, Plots
using Oceananigans
using Oceananigans.Fields
using Oceananigans.AbstractOperations: volume

saved_output_prefix = "horizontal_convection"
saved_output_filename = saved_output_prefix * ".jld2"

# Open the file with our data
file = jldopen(saved_output_filename)

κ = file["closure/κ/b"]

s_timeseries = FieldTimeSeries(saved_output_filename, "s")
b_timeseries = FieldTimeSeries(saved_output_filename, "b")
ζ_timeseries = FieldTimeSeries(saved_output_filename, "ζ")

times = b_timeseries.times

# Coordinate arrays
xs, ys, zs = nodes(s_timeseries[1])
xb, yb, zb = nodes(b_timeseries[1])
xζ, yζ, zζ = nodes(ζ_timeseries[1])
WARNING: using Plots.grid in module Main conflicts with an existing identifier.

Now we're ready to animate.

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

anim = @animate for i in 1:length(times)

    # Load fields from `FieldTimeSeries` and compute χ
    t = times[i]

    s_snapshot = interior(s_timeseries[i], :, 1, :)
    ζ_snapshot = interior(ζ_timeseries[i], :, 1, :)

    b = b_timeseries[i]
    χ = Field(κ * (∂x(b)^2 + ∂z(b)^2))
    compute!(χ)

    b_snapshot = interior(b, :, 1, :)
    χ_snapshot = interior(χ, :, 1, :)

    # determine colorbar limits and contour levels
    slim = 0.6
    slevels = vcat(range(0, stop=slim, length=31), [slim])

    blim = 0.6
    blevels = vcat([-1], range(-blim, stop=blim, length=51), [1])

    ζlim = 9
    ζlevels = range(-ζlim, stop=ζlim, length=41)

    χlim = 0.025
    χlevels = range(0, stop=χlim, length=41)

    @info @sprintf("Drawing frame %d:", i)

    kwargs = (      xlims = (-Lx/2, Lx/2),
                    ylims = (-H, 0),
                   xlabel = "x / H",
                   ylabel = "z / H",
              aspectratio = 1,
                linewidth = 0)

    s_plot = contourf(xs, zs, s_snapshot';
                      clims = (0, slim), levels = slevels,
                      color = :speed, kwargs...)
    s_title = @sprintf("speed √[(u²+w²)/(b⋆H)] @ t=%1.2f", t)

    b_plot = contourf(xb, zb, b_snapshot';
                      clims = (-blim, blim), levels = blevels,
                      color = :thermal, kwargs...)
    b_title = @sprintf("buoyancy, b/b⋆ @ t=%1.2f", t)

    ζ_plot = contourf(xζ, zζ, clamp.(ζ_snapshot', -ζlim, ζlim);
                      clims=(-ζlim, ζlim), levels = ζlevels,
                      color = :balance, kwargs...)
    ζ_title = @sprintf("vorticity, (∂u/∂z - ∂w/∂x) √(H/b⋆) @ t=%1.2f", t)

    χ_plot = contourf(xs, zs, χ_snapshot';
                      clims = (0, χlim), levels = χlevels,
                      color = :dense, kwargs...)
    χ_title = @sprintf("buoyancy dissipation, κ|∇b|² √(H/b⋆⁵) @ t=%1.2f", t)

    plot(s_plot, b_plot, ζ_plot, χ_plot,
           size = (700, 1200),
           link = :x,
         layout = Plots.grid(4, 1),
          title = [s_title b_title ζ_title χ_title])
end

At higher Rayleigh numbers the flow becomes much more vigorous. See, for example, an animation of the voricity of the fluid at $Ra = 10^{12}$ on vimeo.

The Nusselt number

Often we are interested on how much the flow enhances mixing. This is quantified by the Nusselt number, which measures how much the flow enhances mixing compared if only diffusion was in operation. The Nusselt number is given by

\[Nu = \frac{\langle \chi \rangle}{\langle \chi_{\rm diff} \rangle} \, ,\]

where angle brackets above denote both a volume and time average and $\chi_{\rm diff}$ is the buoyancy dissipation that we get without any flow, i.e., the buoyancy dissipation associated with the buoyancy distribution that satisfies

\[\kappa \nabla^2 b_{\rm diff} = 0 \, ,\]

with the same boundary conditions same as our setup. In this case we can readily find that

\[b_{\rm diff}(x, z) = b_s(x) \frac{\cosh \left [2 \pi (H + z) / L_x \right ]}{\cosh(2 \pi H / L_x)} \, ,\]

which implies $\langle \chi_{\rm diff} \rangle = \kappa b_*^2 \pi \tanh(2 \pi Η /Lx)$.

We use the loaded FieldTimeSeries to compute the Nusselt number from buoyancy and the volume average kinetic energy of the fluid.

First we compute the diffusive buoyancy dissipation, $\chi_{\rm diff}$ (which is just a scalar):

H  = file["grid/Lz"]
Lx = file["grid/Lx"]
κ  = file["closure/κ/b"]

χ_diff = κ * b★^2 * π * tanh(2π * H / Lx)

We then create two ReducedFields to store the volume integrals of the kinetic energy density and the buoyancy dissipation. We need the grid to do so; the grid can be recoverd from any FieldTimeSeries, e.g.,

grid = b_timeseries.grid

∫ⱽ_s² = Field{Nothing, Nothing, Nothing}(grid)
∫ⱽ_mod²_∇b = Field{Nothing, Nothing, Nothing}(grid)
1×1×1 Field{Nothing, Nothing, Nothing} reduced over dims = (1, 2, 3) on RectilinearGrid on CPU
├── grid: 128×1×64 RectilinearGrid{Float64, Bounded, Flat, Bounded} on CPU with 3×0×3 halo
├── boundary conditions: FieldBoundaryConditions
│   └── west: Nothing, east: Nothing, south: Nothing, north: Nothing, bottom: Nothing, top: Nothing, immersed: ZeroFlux
└── data: 1×1×1 OffsetArray(::Array{Float64, 3}, 1:1, 1:1, 1:1) with eltype Float64 with indices 1:1×1:1×1:1
    └── max=0.0, min=0.0, mean=0.0

We recover the time from the saved FieldTimeSeries and construct two empty arrays to store the volume-averaged kinetic energy and the instantaneous Nusselt number,

t = b_timeseries.times

kinetic_energy, Nu = zeros(length(t)), zeros(length(t))

Now we can loop over the fields in the FieldTimeSeries, compute KineticEnergy and Nu, and plot.

for i = 1:length(t)
    s = s_timeseries[i]
    sum!(∫ⱽ_s², s^2 * volume)
    kinetic_energy[i] = 0.5 * ∫ⱽ_s²[1, 1, 1]  / (Lx * H)

    b = b_timeseries[i]
    sum!(∫ⱽ_mod²_∇b, (∂x(b)^2 + ∂z(b)^2) * volume)
    Nu[i] = (κ *  ∫ⱽ_mod²_∇b[1, 1, 1]) / χ_diff
end

p1 = plot(t, kinetic_energy,
          xlabel = "time",
          ylabel = "KE / (b⋆H)",
       linewidth = 3,
          legend = :none)

p2 = plot(t, Nu,
          xlabel = "time",
          ylabel = "Nu",
       linewidth = 3,
          legend = :none)

plot(p1, p2, layout = Plots.grid(2, 1))

This page was generated using Literate.jl.