Splitting supercell

This example simulates the development of a splitting supercell thunderstorm, following the idealized test case described by Klemp et al. (2015) and the DCMIP2016 supercell intercomparison by Zarzycki et al. (2019). This benchmark evaluates the model's ability to capture deep moist convection with warm-rain microphysics and strong updrafts.

For microphysics we use the Kessler scheme, which includes prognostic cloud water and rain water with autoconversion, accretion, rain evaporation, and sedimentation processes. This is the same scheme used in the DCMIP2016 supercell intercomparison (Zarzycki et al., 2019).

Physical setup

The simulation initializes a conditionally unstable atmosphere with a warm bubble perturbation that triggers deep convection. The environment includes:

  • A realistic tropospheric potential temperature profile with a tropopause at 12 km
  • Moisture that decreases with height, with relative humidity dropping above the tropopause
  • Wind shear in the lower 5 km to promote storm rotation and supercell development

Potential temperature profile

The background potential temperature follows a piecewise profile (Equation 14 in Klemp et al. (2015)):

\[θ(z) = \begin{cases} θ_0 + (θ_{\rm tr} - θ_0) \left(\frac{z}{z_{\rm tr}}\right)^{5/4} & z \leq z_{\rm tr} \\ θ_{\rm tr} \exp\left(\frac{g}{c_p^d T_{\rm tr}} (z - z_{\rm tr})\right) & z > z_{\rm tr} \end{cases}\]

where $θ_0 = 300 \, {\rm K}$ is the surface potential temperature, $θ_{\rm tr} = 343 \, {\rm K}$ is the tropopause potential temperature, $z_{\rm tr} = 12 \, {\rm km}$ is the tropopause height, and $T_{\rm tr} = 213 \, {\rm K}$ is the tropopause temperature.

Warm bubble perturbation

A localized warm bubble triggers convection (Equations 17–18 in Klemp et al. (2015)):

\[θ'(x, y, z) = \begin{cases} Δθ \cos^2\left(\frac{\pi}{2} R\right) & R < 1 \\ 0 & R \geq 1 \end{cases}\]

where $R = \sqrt{(r/r_h)^2 + ((z-z_c)/r_z)^2}$ is the normalized radius, $r = \sqrt{(x-x_c)^2 + (y-y_c)^2}$ is the horizontal distance from the bubble center, $Δθ = 3 \, {\rm K}$ is the perturbation amplitude, $r_h = 10 \, {\rm km}$ is the horizontal radius, and $r_z = 1.5 \, {\rm km}$ is the vertical radius.

Wind shear profile

The zonal wind increases linearly with height up to the shear layer $z_s = 5 \, {\rm km}$, with a smooth transition zone, providing the environmental shear necessary for supercell development and mesocyclone formation (Equations 15-16 in Klemp et al. (2015)).

using Breeze
using Breeze: DCMIP2016KesslerMicrophysics, TetensFormula
using Oceananigans: Oceananigans
using Oceananigans.Units
using Oceananigans.Grids: znodes

using CairoMakie
using CUDA
using Printf

Domain and grid

The domain is 168 km × 168 km × 20 km with 168 × 168 × 40 grid points, giving 1 km horizontal resolution and 500 m vertical resolution. The grid uses periodic lateral boundary conditions and bounded top/bottom boundaries.

Oceananigans.defaults.FloatType = Float32

Nx, Ny, Nz = 168, 168, 40
Lx, Ly, Lz = 168kilometers, 168kilometers, 20kilometers

grid = RectilinearGrid(GPU(),
                       size = (Nx, Ny, Nz),
                       x = (0, Lx),
                       y = (0, Ly),
                       z = (0, Lz),
                       halo = (5, 5, 5),
                       topology = (Periodic, Periodic, Bounded))
168×168×40 RectilinearGrid{Float32, Periodic, Periodic, Bounded} on CUDAGPU with 5×5×5 halo
├── Periodic x ∈ [0.0, 168000.0) regularly spaced with Δx=1000.0
├── Periodic y ∈ [0.0, 168000.0) regularly spaced with Δy=1000.0
└── Bounded  z ∈ [0.0, 20000.0]  regularly spaced with Δz=500.0

Reference state and dynamics

We define the anelastic reference state with surface pressure $p_0 = 1000 \, {\rm hPa}$ and reference potential temperature $θ_0 = 300 \, {\rm K}$.

constants = ThermodynamicConstants(saturation_vapor_pressure = TetensFormula())

reference_state = ReferenceState(grid, constants,
                                 surface_pressure = 100000,
                                 potential_temperature = 300)

dynamics = AnelasticDynamics(reference_state)
AnelasticDynamics(p₀=100000.0, θ₀=300.0)
└── pressure_anomaly: not materialized

Background atmosphere profiles

The atmospheric stratification parameters define the troposphere-stratosphere transition.

θ₀ = 300       # K - surface potential temperature
θᵖ = 343       # K - tropopause potential temperature
zᵖ = 12000     # m - tropopause height
Tᵖ = 213       # K - tropopause temperature

Wind shear parameters control the low-level environmental wind profile:

zˢ = 5kilometers  # m - shear layer height
uˢ = 30           # m/s - maximum shear wind speed
uᶜ = 15           # m/s - storm motion (Galilean translation speed)

Extract thermodynamic constants for profile calculations:

g = constants.gravitational_acceleration
cᵖᵈ = constants.dry_air.heat_capacity

Background potential temperature profile (Equation 14 in Klemp et al. (2015)):

function θ_background(z)
    θᵗ = θ₀ + (θᵖ - θ₀) * (z / zᵖ)^(5/4)
    θˢ = θᵖ * exp(g / (cᵖᵈ * Tᵖ) * (z - zᵖ))
    return (z <= zᵖ) * θᵗ + (z > zᵖ) * θˢ
end
θ_background (generic function with 1 method)

Relative humidity profile (decreases with height, 25% above tropopause):

ℋ_background(z) = (1 - 3/4 * (z / zᵖ)^(5/4)) * (z <= zᵖ) + 1/4 * (z > zᵖ)
ℋ_background (generic function with 1 method)

Zonal wind profile with linear shear below $zˢ$ and smooth transition (Equations 15-16):

function u_background(z)
    uˡ = uˢ * (z / zˢ) - uᶜ
    uᵗ = (-4/5 + 3 * (z / zˢ) - 5/4 * (z / zˢ)^2) * uˢ - uᶜ
    uᵘ = uˢ - uᶜ
    return (z < (zˢ - 1000)) * uˡ +
           (abs(z - zˢ) <= 1000) * uᵗ +
           (z > (zˢ + 1000)) * uᵘ
end
u_background (generic function with 1 method)

Warm bubble perturbation

The warm bubble parameters following Equations 17–18 in Klemp et al. (2015):

Δθ = 3              # K - perturbation amplitude
rᵇʰ = 10kilometers  # m - bubble horizontal radius
rᵇᵛ = 1500          # m - bubble vertical radius
zᵇ = 1500           # m - bubble center height
xᵇ = Lx / 2         # m - bubble center x-coordinate
yᵇ = Ly / 2         # m - bubble center y-coordinate

The total initial potential temperature combines the background profile with the cosine-squared warm bubble perturbation:

function θᵢ(x, y, z)
    θ̄ = θ_background(z)
    r = sqrt((x - xᵇ)^2 + (y - yᵇ)^2)
    R = sqrt((r / rᵇʰ)^2 + ((z - zᵇ) / rᵇᵛ)^2)
    θ′ = ifelse(R < 1, Δθ * cos((π / 2) * R)^2, 0.0)
    return θ̄ + θ′
end

uᵢ(x, y, z) = u_background(z)
uᵢ (generic function with 1 method)

Visualization of initial conditions and warm bubble perturbation

We visualize the background potential temperature, relative humidity, and wind shear profiles that define the environmental stratification:

θ_profile = set!(Field{Nothing, Nothing, Center}(grid), z -> θ_background(z))
ℋ_profile = set!(Field{Nothing, Nothing, Center}(grid), z -> ℋ_background(z) * 100)
u_profile = set!(Field{Nothing, Nothing, Center}(grid), z -> u_background(z))

fig = Figure(size=(1000, 400), fontsize=14)

axθ = Axis(fig[1, 1], xlabel="θ (K)", ylabel="z (km)", title="Potential temperature")
lines!(axθ, θ_profile, linewidth=2, color=:magenta)
hlines!(axθ, [zᵖ / 1000], color=:gray, linestyle=:dash)

axℋ = Axis(fig[1, 2], xlabel="ℋ (%)", ylabel="z (km)", title="Relative humidity")
lines!(axℋ, ℋ_profile, linewidth=2, color=:dodgerblue)
hlines!(axℋ, [zᵖ / 1000], color=:gray, linestyle=:dash)

axu = Axis(fig[1, 3], xlabel="u (m/s)", ylabel="z (km)", title="Wind profile")
lines!(axu, u_profile, linewidth=2, color=:orangered)
hlines!(axu, [zˢ / 1000], color=:gray, linestyle=:dash)
vlines!(axu, [0], color=:black, linestyle=:dot)

fig

Visualize the warm bubble perturbation on a vertical slice through the domain center:

θ′_slice = set!(Field{Center, Nothing, Center}(grid), (x, z) -> θᵢ(x, yᵇ, z) - θ_background(z))

fig = Figure(size=(700, 400), fontsize=14)
ax = Axis(fig[1, 1], xlabel="x (km)", ylabel="z (km)",
          title="Warm bubble perturbation θ′")

hm = heatmap!(ax, θ′_slice, colormap=:thermal, colorrange=(0, Δθ))
Colorbar(fig[1, 2], hm, label="θ′ (K)")

fig

Model setup

We use the DCMIP2016 Kessler microphysics scheme with high-order WENO advection. The Kessler scheme includes prognostic cloud water and rain water with autoconversion, accretion, rain evaporation, and sedimentation processes.

microphysics = DCMIP2016KesslerMicrophysics()
advection = WENO(order=9, minimum_buffer_upwind_order=3)

model = AtmosphereModel(grid; dynamics, microphysics, advection, thermodynamic_constants=constants)
AtmosphereModel{GPU, RectilinearGrid}(time = 0 seconds, iteration = 0)
├── grid: 168×168×40 RectilinearGrid{Float32, Periodic, Periodic, Bounded} on CUDAGPU with 5×5×5 halo
├── dynamics: AnelasticDynamics(p₀=100000.0, θ₀=300.0)
├── formulation: LiquidIcePotentialTemperatureFormulation
├── thermodynamic_constants: ThermodynamicConstants{Float32}
├── timestepper: SSPRungeKutta3
├── advection scheme: 
│   ├── momentum: WENO{5, Float32, Float32}(order=9)
│   ├── ρθ: WENO{5, Float32, Float32}(order=9)
│   ├── ρqᵗ: WENO{5, Float32, Float32}(order=9)
│   ├── ρqᶜˡ: WENO{5, Float32, Float32}(order=9)
│   └── ρqʳ: WENO{5, Float32, Float32}(order=9)
├── forcing: @NamedTuple{ρu::Returns{Float32}, ρv::Returns{Float32}, ρw::Returns{Float32}, ρθ::Returns{Float32}, ρqᵗ::Returns{Float32}, ρqᶜˡ::Returns{Float32}, ρqʳ::Returns{Float32}, ρe::Returns{Float32}}
├── tracers: ()
├── coriolis: Nothing
└── microphysics: DCMIP2016KesslerMicrophysics

Model initialization

We initialize the model with the previously described initial conditions, including a warm-bubble perturbation. We precompute the RH field to ensure GPU compatibility.

ℋᵢ = set!(CenterField(grid), (x, y, z) -> ℋ_background(z))

set!(model, θ=θᵢ, ℋ=ℋᵢ, u=uᵢ)

Simulation

Run for 2 hours with adaptive time stepping (CFL = 0.7):

simulation = Simulation(model; Δt=2, stop_time=2hours)
conjure_time_step_wizard!(simulation, cfl=0.7)

Output and progress

We set up callbacks to monitor simulation health and collect diagnostics. The maximum vertical velocity is tracked during the simulation to avoid saving large 3D datasets.

θˡⁱ = liquid_ice_potential_temperature(model)
qᶜˡ = model.microphysical_fields.qᶜˡ
qʳ = model.microphysical_fields.qʳ
qᵛ = model.microphysical_fields.qᵛ
u, v, w = model.velocities

wall_clock = Ref(time_ns())

function progress(sim)
    elapsed = 1e-9 * (time_ns() - wall_clock[])

    msg = @sprintf("Iter: %d, t: %s, Δt: %s, wall time: %s, max|u|: %.2f m/s, max w: %.2f m/s, min w: %.2f m/s",
                   iteration(sim), prettytime(sim), prettytime(sim.Δt), prettytime(elapsed),
                   maximum(abs, u), maximum(w), minimum(w))

    msg *= @sprintf(", max(qᵛ): %.2e, max(qᶜˡ): %.2e, max(qʳ): %.2e",
                    maximum(qᵛ), maximum(qᶜˡ), maximum(qʳ))
    @info msg

    return nothing
end

add_callback!(simulation, progress, IterationInterval(100))

Collect maximum vertical velocity time series during simulation:

max_w_ts = []
max_w_times = []

function collect_max_w(sim)
    push!(max_w_times, time(sim))
    push!(max_w_ts, maximum(w))
    return nothing
end

add_callback!(simulation, collect_max_w, TimeInterval(1minutes))

Save horizontal slices at z ≈ 5 km for animation:

z = znodes(grid, Center())
k_5km = searchsortedfirst(z, 5000)
@info "Saving xy slices at z = $(z[k_5km]) m (k = $k_5km)"

slice_outputs = (
    wxy = view(w, :, :, k_5km),
    qʳxy = view(qʳ, :, :, k_5km),
    qᶜˡxy = view(qᶜˡ, :, :, k_5km),
)

slices_filename = "splitting_supercell_slices.jld2"
simulation.output_writers[:slices] = JLD2Writer(model, slice_outputs; filename=slices_filename,
                                                schedule = TimeInterval(2minutes),
                                                overwrite_existing = true)

run!(simulation)
[ Info: Saving xy slices at z = 5250.0 m (k = 11)
[ Info: Initializing simulation...
[ Info: Iter: 0, t: 0 seconds, Δt: 2.200 seconds, wall time: 40.624 seconds, max|u|: 15.00 m/s, max w: 0.00 m/s, min w: -0.00 m/s, max(qᵛ): 2.01e-02, max(qᶜˡ): 0.00e+00, max(qʳ): 0.00e+00
[ Info:     ... simulation initialization complete (54.863 seconds)
[ Info: Executing initial time step...
[ Info:     ... initial time step complete (4.217 seconds).
[ Info: Iter: 100, t: 5.692 minutes, Δt: 5.706 seconds, wall time: 1.235 minutes, max|u|: 18.28 m/s, max w: 6.22 m/s, min w: -0.71 m/s, max(qᵛ): 2.03e-02, max(qᶜˡ): 1.60e-03, max(qʳ): 4.08e-05
[ Info: Iter: 200, t: 15.421 minutes, Δt: 4.049 seconds, wall time: 1.330 minutes, max|u|: 50.26 m/s, max w: 64.92 m/s, min w: -20.73 m/s, max(qᵛ): 2.14e-02, max(qᶜˡ): 6.68e-03, max(qʳ): 1.76e-02
[ Info: Iter: 300, t: 20.413 minutes, Δt: 2.561 seconds, wall time: 1.412 minutes, max|u|: 70.81 m/s, max w: 84.62 m/s, min w: -42.26 m/s, max(qᵛ): 2.12e-02, max(qᶜˡ): 1.59e-02, max(qʳ): 2.25e-02
[ Info: Iter: 400, t: 24.676 minutes, Δt: 2.530 seconds, wall time: 1.495 minutes, max|u|: 133.63 m/s, max w: 119.91 m/s, min w: -47.66 m/s, max(qᵛ): 2.10e-02, max(qᶜˡ): 1.90e-02, max(qʳ): 2.38e-02
[ Info: Iter: 500, t: 28.598 minutes, Δt: 2.383 seconds, wall time: 1.579 minutes, max|u|: 151.41 m/s, max w: 128.51 m/s, min w: -53.50 m/s, max(qᵛ): 2.09e-02, max(qᶜˡ): 1.76e-02, max(qʳ): 2.83e-02
[ Info: Iter: 600, t: 32.327 minutes, Δt: 2.156 seconds, wall time: 1.659 minutes, max|u|: 166.54 m/s, max w: 119.15 m/s, min w: -51.89 m/s, max(qᵛ): 2.10e-02, max(qᶜˡ): 1.16e-02, max(qʳ): 3.14e-02
[ Info: Iter: 700, t: 36.084 minutes, Δt: 2.427 seconds, wall time: 1.746 minutes, max|u|: 175.67 m/s, max w: 106.97 m/s, min w: -55.00 m/s, max(qᵛ): 2.07e-02, max(qᶜˡ): 1.23e-02, max(qʳ): 3.38e-02
[ Info: Iter: 800, t: 40.221 minutes, Δt: 2.694 seconds, wall time: 1.837 minutes, max|u|: 164.78 m/s, max w: 90.24 m/s, min w: -53.22 m/s, max(qᵛ): 2.05e-02, max(qᶜˡ): 1.11e-02, max(qʳ): 3.17e-02
[ Info: Iter: 900, t: 44.722 minutes, Δt: 2.618 seconds, wall time: 1.928 minutes, max|u|: 134.43 m/s, max w: 84.97 m/s, min w: -57.20 m/s, max(qᵛ): 2.05e-02, max(qᶜˡ): 1.29e-02, max(qʳ): 2.78e-02
[ Info: Iter: 1000, t: 49.217 minutes, Δt: 3.119 seconds, wall time: 2.018 minutes, max|u|: 120.43 m/s, max w: 80.73 m/s, min w: -58.45 m/s, max(qᵛ): 2.08e-02, max(qᶜˡ): 1.06e-02, max(qʳ): 2.44e-02
[ Info: Iter: 1100, t: 54.549 minutes, Δt: 3.446 seconds, wall time: 2.105 minutes, max|u|: 90.91 m/s, max w: 70.55 m/s, min w: -51.37 m/s, max(qᵛ): 2.09e-02, max(qᶜˡ): 1.02e-02, max(qʳ): 1.97e-02
[ Info: Iter: 1200, t: 1.011 hours, Δt: 4.149 seconds, wall time: 2.192 minutes, max|u|: 77.34 m/s, max w: 74.37 m/s, min w: -50.93 m/s, max(qᵛ): 2.08e-02, max(qᶜˡ): 6.43e-03, max(qʳ): 1.53e-02
[ Info: Iter: 1300, t: 1.120 hours, Δt: 4.280 seconds, wall time: 2.278 minutes, max|u|: 67.17 m/s, max w: 66.91 m/s, min w: -43.63 m/s, max(qᵛ): 2.04e-02, max(qᶜˡ): 4.72e-03, max(qʳ): 1.69e-02
[ Info: Iter: 1400, t: 1.231 hours, Δt: 4.030 seconds, wall time: 2.366 minutes, max|u|: 65.65 m/s, max w: 71.69 m/s, min w: -52.50 m/s, max(qᵛ): 2.07e-02, max(qᶜˡ): 6.81e-03, max(qʳ): 1.75e-02
[ Info: Iter: 1500, t: 1.343 hours, Δt: 3.882 seconds, wall time: 2.494 minutes, max|u|: 58.91 m/s, max w: 60.26 m/s, min w: -30.39 m/s, max(qᵛ): 2.10e-02, max(qᶜˡ): 4.52e-03, max(qʳ): 1.39e-02
[ Info: Iter: 1600, t: 1.444 hours, Δt: 4.063 seconds, wall time: 2.630 minutes, max|u|: 55.16 m/s, max w: 77.54 m/s, min w: -37.10 m/s, max(qᵛ): 2.14e-02, max(qᶜˡ): 8.55e-03, max(qʳ): 1.56e-02
[ Info: Iter: 1700, t: 1.548 hours, Δt: 3.822 seconds, wall time: 2.769 minutes, max|u|: 65.72 m/s, max w: 67.25 m/s, min w: -56.83 m/s, max(qᵛ): 2.14e-02, max(qᶜˡ): 6.11e-03, max(qʳ): 1.65e-02
[ Info: Iter: 1800, t: 1.657 hours, Δt: 4.007 seconds, wall time: 2.906 minutes, max|u|: 58.04 m/s, max w: 69.50 m/s, min w: -45.03 m/s, max(qᵛ): 2.12e-02, max(qᶜˡ): 7.14e-03, max(qʳ): 1.77e-02
[ Info: Iter: 1900, t: 1.776 hours, Δt: 4.639 seconds, wall time: 3.052 minutes, max|u|: 58.02 m/s, max w: 60.82 m/s, min w: -37.61 m/s, max(qᵛ): 2.12e-02, max(qᶜˡ): 5.76e-03, max(qʳ): 1.68e-02
[ Info: Iter: 2000, t: 1.893 hours, Δt: 4.279 seconds, wall time: 3.196 minutes, max|u|: 63.83 m/s, max w: 72.39 m/s, min w: -53.48 m/s, max(qᵛ): 2.13e-02, max(qᶜˡ): 9.97e-03, max(qʳ): 1.74e-02
[ Info: Iter: 2100, t: 1.998 hours, Δt: 3.948 seconds, wall time: 3.337 minutes, max|u|: 62.46 m/s, max w: 72.88 m/s, min w: -40.58 m/s, max(qᵛ): 2.21e-02, max(qᶜˡ): 6.66e-03, max(qʳ): 1.73e-02
[ Info: Simulation is stopping after running for 3.172 minutes.
[ Info: Simulation time 2 hours equals or exceeds stop time 2 hours.

Animation: horizontal slices at z ≈ 5 km

We create a 3-panel animation showing the storm structure at mid-levels:

  • Vertical velocity $w$: reveals the updraft/downdraft structure
  • Cloud liquid $qᶜˡ$: shows the cloud boundaries
  • Rain $qʳ$: indicates precipitation regions

The simulated supercell exhibits splitting behavior, with the initial storm dividing into right-moving and left-moving cells, consistent with the DCMIP2016 intercomparison results (Zarzycki et al., 2019).

wxy_ts = FieldTimeSeries(slices_filename, "wxy")
qʳxy_ts = FieldTimeSeries(slices_filename, "qʳxy")
qᶜˡxy_ts = FieldTimeSeries(slices_filename, "qᶜˡxy")

times = wxy_ts.times
Nt = length(times)

wlim = maximum(abs, wxy_ts) / 2
qʳlim = maximum(qʳxy_ts) / 4
qᶜˡlim = maximum(qᶜˡxy_ts) / 4

fig = Figure(size=(900, 400), fontsize=12)

axw = Axis(fig[1, 1], aspect=1, xlabel="x (m)", ylabel="y (m)", title="w (m/s)")
axqᶜˡ = Axis(fig[1, 2], aspect=1, xlabel="x (m)", ylabel="y (m)", title="qᶜˡ (kg/kg)")
axqʳ = Axis(fig[1, 3], aspect=1, xlabel="x (m)", ylabel="y (m)", title="qʳ (kg/kg)")

n = Observable(1)
wxy_n = @lift wxy_ts[$n]
qᶜˡxy_n = @lift qᶜˡxy_ts[$n]
qʳxy_n = @lift qʳxy_ts[$n]
title = @lift "Splitting supercell at z ≈ 5 km, t = " * prettytime(times[$n])

hmw = heatmap!(axw, wxy_n, colormap=:balance, colorrange=(-wlim, wlim))
hmqᶜˡ = heatmap!(axqᶜˡ, qᶜˡxy_n, colormap=:dense, colorrange=(0, qᶜˡlim))
hmqʳ = heatmap!(axqʳ, qʳxy_n, colormap=:amp, colorrange=(0, qʳlim))

Colorbar(fig[2, 1], hmw, vertical=false)
Colorbar(fig[2, 2], hmqᶜˡ, vertical=false)
Colorbar(fig[2, 3], hmqʳ, vertical=false)

fig[0, :] = Label(fig, title, fontsize=14, tellwidth=false)

CairoMakie.record(fig, "splitting_supercell_slices.mp4", 1:Nt, framerate=10) do nn
    n[] = nn
end

Results: maximum vertical velocity time series

The maximum updraft velocity is a key diagnostic for supercell intensity. Strong supercells typically develop updrafts exceeding 30–50 m/s.

Our simulated storm intensity is notably stronger than the DCMIP2016 intercomparison results reported by Zarzycki et al. (2019). One explanation is that no explicit numerical diffusion is applied in this simulation. As noted by Klemp et al. (2015), the simulated storm intensity and structure are highly sensitive to numerical diffusion.

fig = Figure(size=(700, 400), fontsize=14)
ax = Axis(fig[1, 1], xlabel="Time (s)", ylabel="Maximum w (m/s)", title="Maximum Vertical Velocity",
          xticks=0:1800:7200)
lines!(ax, max_w_times, max_w_ts, linewidth=2)

fig

Julia version and environment information

This example was executed with the following version of Julia:

using InteractiveUtils: versioninfo
versioninfo()
Julia Version 1.12.4
Commit 01a2eadb047 (2026-01-06 16:56 UTC)
Build Info:
  Official https://julialang.org release
Platform Info:
  OS: Linux (x86_64-linux-gnu)
  CPU: 8 × AMD EPYC 7R13 Processor
  WORD_SIZE: 64
  LLVM: libLLVM-18.1.7 (ORCJIT, znver3)
  GC: Built with stock GC
Threads: 1 default, 1 interactive, 1 GC (on 8 virtual cores)
Environment:
  JULIA_GPG = 3673DF529D9049477F76B37566E3C7DC03D6E495
  JULIA_LOAD_PATH = :@breeze
  JULIA_VERSION = 1.12.4
  JULIA_DEPOT_PATH = /usr/local/share/julia:
  JULIA_PATH = /usr/local/julia
  JULIA_PROJECT = @breeze

These were the top-level packages installed in the environment:

import Pkg
Pkg.status()
Status `/__w/Breeze.jl/Breeze.jl/docs/Project.toml`
  [86bc3604] AtmosphericProfilesLibrary v0.1.7
  [660aa2fb] Breeze v0.3.3 `.`
  [052768ef] CUDA v5.9.6
  [13f3f980] CairoMakie v0.15.8
⌅ [6a9e3e04] CloudMicrophysics v0.29.1
  [e30172f5] Documenter v1.16.1
  [daee34ce] DocumenterCitations v1.4.1
  [98b081ad] Literate v2.21.0
  [85f8d34a] NCDatasets v0.14.11
  [9e8cae18] Oceananigans v0.104.5
  [a01a1ee8] RRTMGP v0.21.7
  [b77e0a4c] InteractiveUtils v1.11.0
  [44cfe95a] Pkg v1.12.1
  [9a3f8284] Random v1.11.0
Info Packages marked with ⌅ have new versions available but compatibility constraints restrict them from upgrading. To see why use `status --outdated`

This page was generated using Literate.jl.