Diurnal cycle of radiative convection

This example simulates a dramatic diurnal cycle of moist convection over a tropical land surface. During the day, the sun heats the ground, driving vigorous boundary-layer convection that lofts moisture into towering cumulus clouds. At night, the surface cools rapidly by longwave emission, stabilizing the boundary layer and shutting off convection.

The key ingredient is a time-varying surface temperature that follows the sun: peaking in the early afternoon and dipping well below the air temperature at night. This creates a strong diurnal contrast — afternoon thunderstorms that die at sunset and don't return until the next morning.

Interactive all-sky RRTMGP radiation computes spectrally-resolved shortwave and longwave fluxes. Saturation-adjustment microphysics diagnoses cloud liquid water that feeds back on the radiation. A stretched vertical grid resolves the cloud layer (100 m spacing below 3 km) while extending to 25 km for a realistic atmospheric column. A stratospheric sponge layer above 8 km prevents spurious temperature drift in the coarse upper cells.

using Breeze
using Oceananigans
using Oceananigans.Units
using Dates: DateTime
using Printf, Random, Statistics
using CairoMakie

using NCDatasets  # Required for RRTMGP lookup tables
using RRTMGP
using CUDA

Random.seed!(2025)
Random.TaskLocalRNG()

Grid

We use a 2D vertical slice (x-z) that is periodic in x and bounded in z. The vertical grid is stretched: fine 100 m cells resolve the cloud layer below 3 km, then a smooth transition to 1 km cells carries the column up to 25 km. This gives RRTMGP a realistic atmospheric column (including the stratosphere) while keeping the total cell count modest.

Nx = 128
Lx = 12800   # 12.8 km

arch = GPU()
Oceananigans.defaults.FloatType = Float32

z = PiecewiseStretchedDiscretization(
    z  = [0, 3000, 8000, 15000],
    Δz = [100,  100, 1000,  1000])

Nz = length(z) - 1

grid = RectilinearGrid(arch;
                       size = (Nx, Nz),
                       x = (0, Lx),
                       z,
                       halo = (5, 5),
                       topology = (Periodic, Flat, Bounded))
128×1×51 RectilinearGrid{Float32, Periodic, Flat, Bounded} on CUDAGPU with 5×0×5 halo
├── Periodic x ∈ [0.0, 12800.0) regularly spaced with Δx=100.0
├── Flat y                      
└── Bounded  z ∈ [0.0, 15000.0] variably spaced with min(Δz)=100.0, max(Δz)=1000.0

Reference state

p₀ = 101325  # Surface pressure [Pa]
θ₀ = 300     # Reference potential temperature [K]

constants = ThermodynamicConstants()

reference_state = ReferenceState(grid, constants;
                                 surface_pressure = p₀,
                                 potential_temperature = θ₀,
                                 vapor_mass_fraction = 0)

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

Background atmosphere

RRTMGP requires trace gas concentrations to compute spectral absorption and emission. We specify well-mixed greenhouse gas concentrations and a tropical ozone profile that transitions from low tropospheric values to a stratospheric peak near 25 km.

@inline function tropical_ozone(z)
    troposphere_O₃ = 30e-9 * (1 + 0.5 * z / 10_000)
    zˢᵗ = 25e3
    Hˢᵗ = 5e3
    stratosphere_O₃ = 8e-6 * exp(-((z - zˢᵗ) / Hˢᵗ)^2)
    χˢᵗ = 1 / (1 + exp(-(z - 15e3) / 2))
    return troposphere_O₃ * (1 - χˢᵗ) + stratosphere_O₃ * χˢᵗ
end

background_atmosphere = BackgroundAtmosphere(
    CO₂ = 348e-6,
    CH₄ = 1650e-9,
    N₂O = 306e-9,
    O₃ = tropical_ozone
)
BackgroundAtmosphere with 6 active gases:
  N₂ = 0.78084
  O₂ = 0.20946
  CO₂ = 348.0 ppm
  CH₄ = 1.65 ppm
  N₂O = 306.0 ppb
  O₃ = tropical_ozone (generic function with 1 method)

Diurnal surface temperature

Over land, the surface temperature swings dramatically with the sun. We model this as a sinusoidal cycle peaking 2 hours after solar noon: Tₛ(t) = T̄ₛ + ΔTₛ cos(2π(t - t_peak) / 24h), with T̄ₛ = 300 K (mean) and ΔTₛ = 10 K (amplitude). This gives 310 K in early afternoon and 290 K at night — a 20 K diurnal range, typical of tropical semi-arid land.

We start at midnight (t = 0), so the surface starts cold (290 K), warms through the morning, peaks at t = 14 h (2 pm local), and cools at night. A Field stores the surface temperature and a callback updates it each time step, keeping both the bulk fluxes and RRTMGP in sync.

T̄ₛ = 300   # Mean surface temperature [K]
ΔTₛ = 20   # Diurnal amplitude [K]

Tₛ = Field{Center, Center, Nothing}(grid)
set!(Tₛ, T̄ₛ - ΔTₛ)  # Start at midnight minimum
128×1×1 Field{Oceananigans.Grids.Center, Oceananigans.Grids.Center, Nothing} reduced over dims = (3,) on Oceananigans.Grids.RectilinearGrid on CUDAGPU
├── grid: 128×1×51 RectilinearGrid{Float32, Periodic, Flat, Bounded} on CUDAGPU with 5×0×5 halo
├── boundary conditions: FieldBoundaryConditions
│   └── west: Periodic, east: Periodic, south: Nothing, north: Nothing, bottom: Nothing, top: Nothing, immersed: Nothing
└── data: 138×1×1 OffsetArray(::CUDA.CuArray{Float32, 3, CUDA.DeviceMemory}, -4:133, 1:1, 1:1) with eltype Float32 with indices -4:133×1:1×1:1
    └── max=280.0, min=280.0, mean=280.0

Radiation with a diurnal cycle

We place the domain at 15°N latitude on the prime meridian and start at midnight on the spring equinox (March 20). The sun rises at t ≈ 6 h, reaches noon at t ≈ 12 h, and sets at t ≈ 18 h.

latitude = 15 # 15°N

radiation = RadiativeTransferModel(grid, AllSkyOptics(), constants;
                                   surface_temperature = Tₛ,
                                   surface_albedo = 0.20,
                                   surface_emissivity = 0.95,
                                   solar_constant = 1361,
                                   background_atmosphere,
                                   coordinate = (0, latitude),
                                   epoch = DateTime(2020, 3, 20, 0, 0, 0),
                                   schedule = TimeInterval(5minutes),
                                   liquid_effective_radius = ConstantRadiusParticles(10e-6),
                                   ice_effective_radius = ConstantRadiusParticles(30e-6))
RadiativeTransferModel
├── solar_constant: 1361.0 W m⁻²
├── surface_temperature: 128×1×1 Field{Oceananigans.Grids.Center, Oceananigans.Grids.Center, Nothing} reduced over dims = (3,) on Oceananigans.Grids.RectilinearGrid on CUDAGPU
├── grid: 128×1×51 RectilinearGrid{Float32, Periodic, Flat, Bounded} on CUDAGPU with 5×0×5 halo
├── boundary conditions: FieldBoundaryConditions
│   └── west: Periodic, east: Periodic, south: Nothing, north: Nothing, bottom: Nothing, top: Nothing, immersed: Nothing
└── data: 138×1×1 OffsetArray(::CUDA.CuArray{Float32, 3, CUDA.DeviceMemory}, -4:133, 1:1, 1:1) with eltype Float32 with indices -4:133×1:1×1:1
    └── max=280.0, min=280.0, mean=280.0 K
├── surface_emissivity: ConstantField(0.95)
├── direct_surface_albedo: ConstantField(0.2)
├── liquid_effective_radius: Breeze.AtmosphereModels.ConstantRadiusParticles{Float32}(1.0f-5)
├── ice_effective_radius: Breeze.AtmosphereModels.ConstantRadiusParticles{Float32}(3.0f-5)
└── diffuse_surface_albedo: ConstantField(0.2)

Surface fluxes

Bulk aerodynamic formulae provide surface sensible heat, moisture, and momentum fluxes driven by the time-varying surface temperature. During the day Tₛ > Tair drives strong upward fluxes; at night Tₛ < Tair can produce downward fluxes that cool the boundary layer.

Cᴰ = Cᵀ = 1e-3
Cᵛ = 1.2e-3
Uᵍ = 1  # Gustiness [m/s]

ρθ_flux = BulkSensibleHeatFlux(coefficient=Cᵀ, gustiness=Uᵍ, surface_temperature=Tₛ)
ρqᵗ_flux = BulkVaporFlux(coefficient=Cᵛ, gustiness=Uᵍ, surface_temperature=Tₛ)

ρθ_bcs = FieldBoundaryConditions(bottom=ρθ_flux)
ρqᵗ_bcs = FieldBoundaryConditions(bottom=ρqᵗ_flux)
ρu_bcs = FieldBoundaryConditions(bottom=Breeze.BulkDrag(coefficient=Cᴰ, gustiness=Uᵍ))
Oceananigans.FieldBoundaryConditions, with boundary conditions
├── west: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── east: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── south: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── north: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── bottom: FluxBoundaryCondition: BulkDragFunction(direction=Nothing, coefficient=0.001, gustiness=1)
├── top: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
└── immersed: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)

Microphysics

Warm-phase saturation adjustment instantly converts supersaturated moisture to cloud liquid water.

microphysics = SaturationAdjustment(equilibrium=WarmPhaseEquilibrium())
Breeze.Microphysics.SaturationAdjustment{Breeze.Thermodynamics.WarmPhaseEquilibrium, Float32}(0.001f0, Inf32, Breeze.Thermodynamics.WarmPhaseEquilibrium())

Stratospheric sponge

The domain extends to 25 km, but the initial stratosphere isn't in radiative equilibrium: ozone absorbs shortwave radiation and the coarse upper cells respond strongly. A Newtonian relaxation of temperature toward the initial profile above 8 km keeps the stratosphere anchored without affecting the tropospheric dynamics. We apply this as an energy forcing on ρe, which Breeze automatically converts to a ρθ tendency.

Tᵣ = reference_state.temperature
ρᵣ = reference_state.density
cᵖᵈ = constants.dry_air.heat_capacity / constants.dry_air.molar_mass  # J/(kg·K)
τ_sponge = 6hours

@inline function stratospheric_relaxation(i, j, k, grid, clock, model_fields, p)
    @inbounds T = model_fields.T[i, j, k]
    @inbounds Tᵣ = p.Tᵣ[i, j, k]
    @inbounds ρ = p.ρᵣ[i, j, k]
    z = znode(i, j, k, grid, Center(), Center(), Center())
    α = clamp((z - 8000) / 4000, 0, 1)
    ∂T∂t = -α * (T - Tᵣ) / p.τ
    return ρ * p.cᵖᵈ * ∂T∂t
end

sponge = Forcing(stratospheric_relaxation; discrete_form=true,
                 parameters=(; Tᵣ, ρᵣ, cᵖᵈ, τ=τ_sponge))

forcing = (; ρe=sponge)
(ρe = DiscreteForcing{@NamedTuple{Tᵣ::Oceananigans.Fields.Field{Nothing, Nothing, Oceananigans.Grids.Center, Nothing, Oceananigans.Grids.RectilinearGrid{Float32, Oceananigans.Grids.Periodic, Oceananigans.Grids.Flat, Oceananigans.Grids.Bounded, Oceananigans.Grids.StaticVerticalDiscretization{OffsetArrays.OffsetVector{Float32, CUDA.CuArray{Float32, 1, CUDA.DeviceMemory}}, OffsetArrays.OffsetVector{Float32, CUDA.CuArray{Float32, 1, CUDA.DeviceMemory}}, OffsetArrays.OffsetVector{Float32, CUDA.CuArray{Float32, 1, CUDA.DeviceMemory}}, OffsetArrays.OffsetVector{Float32, CUDA.CuArray{Float32, 1, CUDA.DeviceMemory}}}, Float32, Float32, OffsetArrays.OffsetVector{Float32, StepRangeLen{Float32, Float64, Float64, Int64}}, Nothing, Oceananigans.Architectures.GPU{CUDA.CUDAKernels.CUDABackend}}, Tuple{Colon, Colon, Colon}, OffsetArrays.OffsetArray{Float32, 3, CUDA.CuArray{Float32, 3, CUDA.DeviceMemory}}, Float32, Oceananigans.BoundaryConditions.FieldBoundaryConditions{Nothing, Nothing, Nothing, Nothing, Oceananigans.BoundaryConditions.NoFluxBoundaryCondition, Oceananigans.BoundaryConditions.NoFluxBoundaryCondition, Nothing, @NamedTuple{bottom_and_top::KernelAbstractions.Kernel{CUDA.CUDAKernels.CUDABackend, KernelAbstractions.NDIteration.StaticSize{(1, 1)}, Oceananigans.Utils.OffsetStaticSize{(1:1, 1:1)}, typeof(Oceananigans.BoundaryConditions.gpu__fill_bottom_and_top_halo!)}, south_and_north::Nothing, west_and_east::Nothing}, @NamedTuple{bottom_and_top::Tuple{Oceananigans.BoundaryConditions.NoFluxBoundaryCondition, Oceananigans.BoundaryConditions.NoFluxBoundaryCondition}, south_and_north::Tuple{Nothing, Nothing}, west_and_east::Tuple{Nothing, Nothing}}}, Nothing, Nothing}, ρᵣ::Oceananigans.Fields.Field{Nothing, Nothing, Oceananigans.Grids.Center, Nothing, Oceananigans.Grids.RectilinearGrid{Float32, Oceananigans.Grids.Periodic, Oceananigans.Grids.Flat, Oceananigans.Grids.Bounded, Oceananigans.Grids.StaticVerticalDiscretization{OffsetArrays.OffsetVector{Float32, CUDA.CuArray{Float32, 1, CUDA.DeviceMemory}}, OffsetArrays.OffsetVector{Float32, CUDA.CuArray{Float32, 1, CUDA.DeviceMemory}}, OffsetArrays.OffsetVector{Float32, CUDA.CuArray{Float32, 1, CUDA.DeviceMemory}}, OffsetArrays.OffsetVector{Float32, CUDA.CuArray{Float32, 1, CUDA.DeviceMemory}}}, Float32, Float32, OffsetArrays.OffsetVector{Float32, StepRangeLen{Float32, Float64, Float64, Int64}}, Nothing, Oceananigans.Architectures.GPU{CUDA.CUDAKernels.CUDABackend}}, Tuple{Colon, Colon, Colon}, OffsetArrays.OffsetArray{Float32, 3, CUDA.CuArray{Float32, 3, CUDA.DeviceMemory}}, Float32, Oceananigans.BoundaryConditions.FieldBoundaryConditions{Nothing, Nothing, Nothing, Nothing, Oceananigans.BoundaryConditions.BoundaryCondition{Oceananigans.BoundaryConditions.Value, Float32}, Oceananigans.BoundaryConditions.NoFluxBoundaryCondition, Nothing, @NamedTuple{bottom_and_top::KernelAbstractions.Kernel{CUDA.CUDAKernels.CUDABackend, KernelAbstractions.NDIteration.StaticSize{(1, 1)}, Oceananigans.Utils.OffsetStaticSize{(1:1, 1:1)}, typeof(Oceananigans.BoundaryConditions.gpu__fill_bottom_and_top_halo!)}, south_and_north::Nothing, west_and_east::Nothing}, @NamedTuple{bottom_and_top::Tuple{Oceananigans.BoundaryConditions.BoundaryCondition{Oceananigans.BoundaryConditions.Value, Float32}, Oceananigans.BoundaryConditions.NoFluxBoundaryCondition}, south_and_north::Tuple{Nothing, Nothing}, west_and_east::Tuple{Nothing, Nothing}}}, Nothing, Nothing}, cᵖᵈ::Float32, τ::Float64}}
├── func: stratospheric_relaxation (generic function with 1 method)
└── parameters: (Tᵣ = 1×1×51 Field{Nothing, Nothing, Oceananigans.Grids.Center} reduced over dims = (1, 2) on Oceananigans.Grids.RectilinearGrid on CUDAGPU
├── grid: 128×1×51 RectilinearGrid{Float32, Periodic, Flat, Bounded} on CUDAGPU with 5×0×5 halo
├── boundary conditions: FieldBoundaryConditions
│   └── west: Nothing, east: Nothing, south: Nothing, north: Nothing, bottom: ZeroFlux, top: ZeroFlux, immersed: Nothing
└── data: 1×1×61 OffsetArray(::CUDA.CuArray{Float32, 3, CUDA.DeviceMemory}, 1:1, 1:1, -4:56) with eltype Float32 with indices 1:1×1:1×-4:56
    └── max=300.64, min=159.059, mean=264.54, ρᵣ = 1×1×51 Field{Nothing, Nothing, Oceananigans.Grids.Center} reduced over dims = (1, 2) on Oceananigans.Grids.RectilinearGrid on CUDAGPU
├── grid: 128×1×51 RectilinearGrid{Float32, Periodic, Flat, Bounded} on CUDAGPU with 5×0×5 halo
├── boundary conditions: FieldBoundaryConditions
│   └── west: Nothing, east: Nothing, south: Nothing, north: Nothing, bottom: Value, top: ZeroFlux, immersed: Nothing
└── data: 1×1×61 OffsetArray(::CUDA.CuArray{Float32, 3, CUDA.DeviceMemory}, 1:1, 1:1, -4:56) with eltype Float32 with indices 1:1×1:1×-4:56
    └── max=1.16764, min=0.237475, mean=0.875196, cᵖᵈ = 34691.06f0, τ = 21600.0),)

Model assembly

coriolis = FPlane(; latitude)
boundary_conditions = (ρθ=ρθ_bcs, ρqᵗ=ρqᵗ_bcs, ρu=ρu_bcs)

weno_order = 5
momentum_advection = WENO(order=weno_order)

scalar_advection = (ρθ  = WENO(order=weno_order),
                    ρqᵗ = WENO(order=weno_order, bounds=(0, 1)))

model = AtmosphereModel(grid; dynamics, microphysics, radiation, forcing,
                        momentum_advection, scalar_advection,
                        boundary_conditions, coriolis)
AtmosphereModel{GPU, RectilinearGrid}(time = 0 seconds, iteration = 0)
├── grid: 128×1×51 RectilinearGrid{Float32, Periodic, Flat, Bounded} on CUDAGPU with 5×0×5 halo
├── dynamics: AnelasticDynamics(p₀=101325.0, θ₀=300.0)
├── formulation: LiquidIcePotentialTemperatureFormulation
├── thermodynamic_constants: ThermodynamicConstants{Float32}
├── timestepper: SSPRungeKutta3
├── advection scheme: 
│   ├── momentum: WENO{3, Float32, Float32}(order=5)
│   ├── ρθ: WENO{3, Float32, Float32}(order=5)
│   └── ρqᵗ: WENO{3, Float32, Float32}(order=5, bounds=(0.0f0, 1.0f0))
├── forcing: ρe=>DiscreteForcing
├── tracers: ()
├── coriolis: FPlane{Float32}(f=3.77468e-5)
└── microphysics: SaturationAdjustment

Initial conditions

The sounding has a dry-adiabatic sub-cloud layer (0–1 km) capped by a conditionally unstable troposphere (5 K/km lapse rate) that transitions to an isothermal stratosphere at 210 K. Moisture is 20 g/kg at the surface with a 2.5 km scale height, typical of the tropical maritime boundary layer.

function Tᵇᵍ(z)
    T₀ = 300
    Tˢᵗ = 210 # stratosphere temperature
    T = T₀ - 1e-3 * max(z, 1000) - 5e-3 * max(0, z - 1000)
    return max(T, Tˢᵗ)
end

uᵢ(x, z) = -5 * max(1 - z / 3000, 0)
uᵢ (generic function with 1 method)

Random perturbations in the lowest 1 km trigger convection.

δT = 2
δℋ = 1e-2
zδ = 1000

ϵ() = rand() - 0.5
Tᵢ(x, z) = Tᵇᵍ(z) + δT * ϵ() * (z < zδ)
ℋᵢ(x, z) = (0.5 + δℋ * ϵ()) * (z < zδ)
ℋᵢ (generic function with 1 method)

After setting initial conditions, we recompute the reference state from the horizontally-averaged model state. This is important for Float32 accuracy in tall domains: the default dry-adiabat reference state diverges from the actual stratospheric profile, causing large density errors that overwhelm Float32 precision. set_to_mean! adjusts ρᵣ to match the current state while rescaling density-weighted prognostic fields to preserve specific quantities.

set!(model; T=Tᵢ, ℋ=ℋᵢ, u=uᵢ)

reference_state = model.dynamics.reference_state
set_to_mean!(reference_state, model, rescale_densities=true)

T = model.temperature
qᵗ = model.specific_moisture
u, w = model.velocities.u, model.velocities.w
qˡ = model.microphysical_fields.qˡ

@info "Diurnal Radiative Convection (2D)"
@info "Grid: $(Nx) × $(Nz) (stretched), domain: $(Lx/1000) km × 25 km"
@info "Initial T range: $(minimum(T)) – $(maximum(T)) K"
@info "Initial qᵗ range: $(minimum(qᵗ)*1000) – $(maximum(qᵗ)*1000) g/kg"
[ Info: Diurnal Radiative Convection (2D)
[ Info: Grid: 128 × 51 (stretched), domain: 12.8 km × 25 km
[ Info: Initial T range: 238.96986 – 300.16586 K
[ Info: Initial qᵗ range: 0.0 – 12.126032 g/kg

Simulation

We run for two full diurnal cycles (48 hours) starting at midnight, so the on/off pattern of convection repeats convincingly.

simulation = Simulation(model; Δt=1, stop_time=3days)
conjure_time_step_wizard!(simulation, cfl=0.7)

Surface temperature callback

At each time step we update the surface temperature field following a cosine curve that peaks 2 hours after solar noon (t = 14 h local). The period is 24 hours and the simulation starts at midnight.

function update_surface_temperature!(sim)
    t = time(sim)
    t_peak = 14hour  # peak at 14:00 local (2 pm)
    Tₛ_now = T̄ₛ + ΔTₛ * cos(2π * (t - t_peak) / day)
    set!(Tₛ, Tₛ_now)
    return nothing
end

add_callback!(simulation, update_surface_temperature!, TimeInterval(1minute))

wall_clock = Ref(time_ns())

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

    wmax = maximum(abs, w)
    Tmin, Tmax = extrema(T)
    qˡmax = maximum(qˡ)
    Tₛ_now = T̄ₛ + ΔTₛ * cos(2π * (time(sim) - 14hours) / 24hours)

    OLR = mean(view(radiation.upwelling_longwave_flux, :, 1, Nz+1))

    msg = @sprintf("Iter: %5d, t: %8s, Δt: %5.1fs, wall: %8s",
                   iteration(sim), prettytime(sim), sim.Δt, prettytime(elapsed))
    msg *= @sprintf(", max|w|: %5.2f m/s, T: [%5.1f, %5.1f] K, max(qˡ): %.2e",
                   wmax, Tmin, Tmax, qˡmax)
    msg *= @sprintf(", Tₛ: %.1f K, OLR: %.1f W/m²", Tₛ_now, OLR)
    @info msg

    wall_clock[] = time_ns()
    return nothing
end

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

Output

Horizontally-averaged profiles are saved every hour (time-averaged) and 2D slices every 10 minutes for animation.

qᵛ = model.microphysical_fields.qᵛ
Fᴿ = radiation.flux_divergence

outputs = (; u, w, T, qˡ, qᵛ, Fᴿ)
avg_outputs = NamedTuple(name => Average(outputs[name], dims=1) for name in keys(outputs))

filename = "radiative_convection"
averages_filename = filename * "_averages.jld2"
slices_filename = filename * "_slices.jld2"

simulation.output_writers[:averages] = JLD2Writer(model, avg_outputs;
                                                  filename = averages_filename,
                                                  schedule = AveragedTimeInterval(1hour),
                                                  overwrite_existing = true)

slice_outputs = (; w, qᵛ, T)
simulation.output_writers[:slices] = JLD2Writer(model, slice_outputs;
                                                filename = slices_filename,
                                                schedule = TimeInterval(10minutes),
                                                overwrite_existing = true)

@info "Starting simulation..."
run!(simulation)
@info "Simulation completed!"
[ Info: Starting simulation...
[ Info: Initializing simulation...
[ Info: Iter:     0, t: 0 seconds, Δt:   1.1s, wall: 3.032 minutes, max|w|:  0.00 m/s, T: [239.0, 300.2] K, max(qˡ): 0.00e+00, Tₛ: 282.7 K, OLR: 374.3 W/m²
[ Info:     ... simulation initialization complete (31.042 seconds)
[ Info: Executing initial time step...
[ Info:     ... initial time step complete (7.545 seconds).
[ Info: Iter:  1000, t: 1.277 hours, Δt:   3.2s, wall: 34.511 seconds, max|w|: 12.91 m/s, T: [217.9, 302.9] K, max(qˡ): 0.00e+00, Tₛ: 280.4 K, OLR: 349.9 W/m²
[ Info: Iter:  2000, t: 2.304 hours, Δt:   4.3s, wall: 18.664 seconds, max|w|: 10.80 m/s, T: [217.5, 302.3] K, max(qˡ): 0.00e+00, Tₛ: 280.1 K, OLR: 341.7 W/m²
[ Info: Iter:  3000, t: 3.433 hours, Δt:   4.6s, wall: 19.434 seconds, max|w|:  7.01 m/s, T: [217.5, 303.6] K, max(qˡ): 0.00e+00, Tₛ: 281.4 K, OLR: 336.4 W/m²
[ Info: Iter:  4000, t: 4.557 hours, Δt:   4.8s, wall: 19.406 seconds, max|w|:  9.35 m/s, T: [217.2, 302.8] K, max(qˡ): 0.00e+00, Tₛ: 284.3 K, OLR: 336.7 W/m²
[ Info: Iter:  5000, t: 5.615 hours, Δt:   3.1s, wall: 18.875 seconds, max|w|: 14.43 m/s, T: [217.3, 304.5] K, max(qˡ): 0.00e+00, Tₛ: 288.3 K, OLR: 339.1 W/m²
[ Info: Iter:  6000, t: 6.634 hours, Δt:   3.8s, wall: 18.558 seconds, max|w|: 11.11 m/s, T: [217.5, 304.8] K, max(qˡ): 0.00e+00, Tₛ: 293.0 K, OLR: 343.3 W/m²
[ Info: Iter:  7000, t: 7.554 hours, Δt:   3.3s, wall: 18.025 seconds, max|w|:  9.85 m/s, T: [217.5, 306.5] K, max(qˡ): 0.00e+00, Tₛ: 297.7 K, OLR: 348.8 W/m²
[ Info: Iter:  8000, t: 8.568 hours, Δt:   3.7s, wall: 18.459 seconds, max|w|: 10.73 m/s, T: [217.4, 303.8] K, max(qˡ): 0.00e+00, Tₛ: 303.0 K, OLR: 356.5 W/m²
[ Info: Iter:  9000, t: 9.636 hours, Δt:   3.9s, wall: 18.808 seconds, max|w|: 11.88 m/s, T: [217.5, 303.9] K, max(qˡ): 0.00e+00, Tₛ: 308.3 K, OLR: 364.7 W/m²
[ Info: Iter: 10000, t: 10.696 hours, Δt:   4.1s, wall: 18.841 seconds, max|w|: 11.01 m/s, T: [217.7, 304.7] K, max(qˡ): 0.00e+00, Tₛ: 313.0 K, OLR: 372.3 W/m²
[ Info: Iter: 11000, t: 11.810 hours, Δt:   4.0s, wall: 18.852 seconds, max|w|: 11.34 m/s, T: [217.4, 304.0] K, max(qˡ): 0.00e+00, Tₛ: 316.8 K, OLR: 377.8 W/m²
[ Info: Iter: 12000, t: 12.950 hours, Δt:   4.0s, wall: 19.160 seconds, max|w|: 12.22 m/s, T: [217.4, 304.4] K, max(qˡ): 0.00e+00, Tₛ: 319.2 K, OLR: 378.7 W/m²
[ Info: Iter: 13000, t: 13.994 hours, Δt:   4.6s, wall: 18.425 seconds, max|w|:  9.31 m/s, T: [217.5, 305.7] K, max(qˡ): 0.00e+00, Tₛ: 320.0 K, OLR: 378.6 W/m²
[ Info: Iter: 14000, t: 15.058 hours, Δt:   5.1s, wall: 18.853 seconds, max|w|:  7.57 m/s, T: [217.6, 305.8] K, max(qˡ): 0.00e+00, Tₛ: 319.2 K, OLR: 375.2 W/m²
[ Info: Iter: 15000, t: 16.140 hours, Δt:   3.4s, wall: 18.782 seconds, max|w|: 13.78 m/s, T: [217.3, 304.3] K, max(qˡ): 1.80e-04, Tₛ: 316.9 K, OLR: 369.1 W/m²
[ Info: Iter: 16000, t: 17.076 hours, Δt:   2.8s, wall: 18.049 seconds, max|w|: 17.17 m/s, T: [217.5, 304.0] K, max(qˡ): 1.12e-03, Tₛ: 313.9 K, OLR: 360.6 W/m²
[ Info: Iter: 17000, t: 18.055 hours, Δt:   3.3s, wall: 18.493 seconds, max|w|: 15.21 m/s, T: [217.5, 303.1] K, max(qˡ): 1.77e-03, Tₛ: 309.8 K, OLR: 352.5 W/m²
[ Info: Iter: 18000, t: 19.118 hours, Δt:   3.0s, wall: 18.836 seconds, max|w|: 12.59 m/s, T: [217.5, 303.3] K, max(qˡ): 3.62e-04, Tₛ: 304.6 K, OLR: 343.1 W/m²
[ Info: Iter: 19000, t: 20.133 hours, Δt:   3.4s, wall: 18.483 seconds, max|w|: 13.32 m/s, T: [217.6, 304.2] K, max(qˡ): 0.00e+00, Tₛ: 299.3 K, OLR: 334.2 W/m²
[ Info: Iter: 20000, t: 21.236 hours, Δt:   4.0s, wall: 18.836 seconds, max|w|: 14.20 m/s, T: [217.3, 304.4] K, max(qˡ): 0.00e+00, Tₛ: 293.6 K, OLR: 324.3 W/m²
[ Info: Iter: 21000, t: 22.111 hours, Δt:   4.1s, wall: 18.039 seconds, max|w|:  9.17 m/s, T: [217.6, 304.0] K, max(qˡ): 0.00e+00, Tₛ: 289.5 K, OLR: 317.4 W/m²
[ Info: Iter: 22000, t: 23.172 hours, Δt:   3.4s, wall: 18.810 seconds, max|w|: 13.93 m/s, T: [217.2, 304.8] K, max(qˡ): 9.08e-04, Tₛ: 285.2 K, OLR: 310.1 W/m²
[ Info: Iter: 23000, t: 1.008 days, Δt:   4.0s, wall: 18.422 seconds, max|w|: 11.67 m/s, T: [217.7, 303.0] K, max(qˡ): 0.00e+00, Tₛ: 282.2 K, OLR: 305.2 W/m²
[ Info: Iter: 24000, t: 1.052 days, Δt:   3.5s, wall: 18.507 seconds, max|w|: 10.46 m/s, T: [217.6, 305.2] K, max(qˡ): 1.69e-04, Tₛ: 280.4 K, OLR: 302.2 W/m²
[ Info: Iter: 25000, t: 1.096 days, Δt:   3.8s, wall: 18.818 seconds, max|w|: 11.83 m/s, T: [217.6, 303.3] K, max(qˡ): 0.00e+00, Tₛ: 280.1 K, OLR: 300.8 W/m²
[ Info: Iter: 26000, t: 1.141 days, Δt:   3.8s, wall: 18.818 seconds, max|w|: 11.20 m/s, T: [217.4, 303.9] K, max(qˡ): 5.12e-04, Tₛ: 281.3 K, OLR: 301.1 W/m²
[ Info: Iter: 27000, t: 1.181 days, Δt:   3.6s, wall: 18.448 seconds, max|w|: 12.64 m/s, T: [217.6, 304.3] K, max(qˡ): 1.88e-03, Tₛ: 283.6 K, OLR: 303.0 W/m²
[ Info: Iter: 28000, t: 1.219 days, Δt:   3.2s, wall: 17.683 seconds, max|w|: 12.24 m/s, T: [217.5, 302.9] K, max(qˡ): 1.51e-03, Tₛ: 286.8 K, OLR: 306.1 W/m²
[ Info: Iter: 29000, t: 1.263 days, Δt:   4.8s, wall: 18.874 seconds, max|w|:  8.31 m/s, T: [217.4, 302.7] K, max(qˡ): 0.00e+00, Tₛ: 291.5 K, OLR: 311.9 W/m²
[ Info: Iter: 30000, t: 1.307 days, Δt:   3.4s, wall: 18.866 seconds, max|w|: 13.48 m/s, T: [217.6, 304.8] K, max(qˡ): 0.00e+00, Tₛ: 296.7 K, OLR: 319.2 W/m²
[ Info: Iter: 31000, t: 1.349 days, Δt:   4.4s, wall: 18.454 seconds, max|w|: 10.16 m/s, T: [217.7, 304.5] K, max(qˡ): 3.48e-04, Tₛ: 302.0 K, OLR: 326.3 W/m²
[ Info: Iter: 32000, t: 1.390 days, Δt:   3.7s, wall: 18.481 seconds, max|w|: 11.22 m/s, T: [217.6, 303.5] K, max(qˡ): 0.00e+00, Tₛ: 307.0 K, OLR: 333.6 W/m²
[ Info: Iter: 33000, t: 1.435 days, Δt:   3.6s, wall: 18.843 seconds, max|w|: 12.05 m/s, T: [217.5, 304.4] K, max(qˡ): 0.00e+00, Tₛ: 312.0 K, OLR: 340.8 W/m²
[ Info: Iter: 34000, t: 1.479 days, Δt:   3.8s, wall: 18.847 seconds, max|w|: 11.38 m/s, T: [217.5, 304.9] K, max(qˡ): 4.92e-04, Tₛ: 315.9 K, OLR: 345.7 W/m²
[ Info: Iter: 35000, t: 1.520 days, Δt:   3.1s, wall: 18.149 seconds, max|w|: 15.06 m/s, T: [217.4, 305.4] K, max(qˡ): 7.06e-04, Tₛ: 318.4 K, OLR: 348.2 W/m²
[ Info: Iter: 36000, t: 1.559 days, Δt:   4.1s, wall: 18.531 seconds, max|w|: 12.24 m/s, T: [217.7, 305.7] K, max(qˡ): 0.00e+00, Tₛ: 319.8 K, OLR: 350.0 W/m²
[ Info: Iter: 37000, t: 1.603 days, Δt:   3.9s, wall: 18.526 seconds, max|w|: 10.44 m/s, T: [217.4, 306.5] K, max(qˡ): 0.00e+00, Tₛ: 319.9 K, OLR: 348.7 W/m²
[ Info: Iter: 38000, t: 1.648 days, Δt:   4.1s, wall: 18.841 seconds, max|w|: 10.72 m/s, T: [217.6, 303.9] K, max(qˡ): 2.16e-03, Tₛ: 318.4 K, OLR: 343.7 W/m²
[ Info: Iter: 39000, t: 1.691 days, Δt:   2.9s, wall: 18.470 seconds, max|w|: 14.41 m/s, T: [217.5, 303.4] K, max(qˡ): 1.61e-03, Tₛ: 315.6 K, OLR: 338.4 W/m²
[ Info: Iter: 40000, t: 1.728 days, Δt:   3.1s, wall: 18.047 seconds, max|w|: 13.94 m/s, T: [217.5, 304.2] K, max(qˡ): 1.77e-03, Tₛ: 312.2 K, OLR: 333.3 W/m²
[ Info: Iter: 41000, t: 1.767 days, Δt:   3.7s, wall: 18.040 seconds, max|w|:  8.31 m/s, T: [217.4, 305.0] K, max(qˡ): 8.57e-04, Tₛ: 308.1 K, OLR: 327.3 W/m²
[ Info: Iter: 42000, t: 1.810 days, Δt:   4.0s, wall: 18.779 seconds, max|w|: 12.81 m/s, T: [217.4, 302.7] K, max(qˡ): 1.27e-03, Tₛ: 302.9 K, OLR: 318.6 W/m²
[ Info: Iter: 43000, t: 1.855 days, Δt:   4.1s, wall: 18.869 seconds, max|w|:  9.94 m/s, T: [217.4, 302.7] K, max(qˡ): 7.40e-04, Tₛ: 297.3 K, OLR: 312.0 W/m²
[ Info: Iter: 44000, t: 1.898 days, Δt:   3.5s, wall: 18.435 seconds, max|w|: 13.17 m/s, T: [217.5, 303.8] K, max(qˡ): 1.03e-03, Tₛ: 292.0 K, OLR: 306.1 W/m²
[ Info: Iter: 45000, t: 1.941 days, Δt:   4.4s, wall: 18.773 seconds, max|w|: 10.67 m/s, T: [217.3, 304.1] K, max(qˡ): 6.06e-04, Tₛ: 287.5 K, OLR: 300.2 W/m²
[ Info: Iter: 46000, t: 1.980 days, Δt:   2.7s, wall: 18.075 seconds, max|w|: 15.48 m/s, T: [217.6, 303.8] K, max(qˡ): 3.90e-03, Tₛ: 284.0 K, OLR: 295.9 W/m²
[ Info: Iter: 47000, t: 2.019 days, Δt:   3.7s, wall: 18.044 seconds, max|w|:  9.53 m/s, T: [217.6, 302.6] K, max(qˡ): 1.39e-03, Tₛ: 281.6 K, OLR: 293.4 W/m²
[ Info: Iter: 48000, t: 2.058 days, Δt:   3.8s, wall: 18.059 seconds, max|w|:  8.57 m/s, T: [217.5, 303.3] K, max(qˡ): 5.89e-04, Tₛ: 280.3 K, OLR: 291.6 W/m²
[ Info: Iter: 49000, t: 2.102 days, Δt:   4.8s, wall: 18.837 seconds, max|w|: 10.01 m/s, T: [217.5, 303.2] K, max(qˡ): 2.54e-03, Tₛ: 280.1 K, OLR: 290.5 W/m²
[ Info: Iter: 50000, t: 2.150 days, Δt:   4.2s, wall: 19.222 seconds, max|w|:  9.86 m/s, T: [217.3, 302.4] K, max(qˡ): 0.00e+00, Tₛ: 281.7 K, OLR: 291.5 W/m²
[ Info: Iter: 51000, t: 2.197 days, Δt:   3.3s, wall: 18.894 seconds, max|w|: 11.44 m/s, T: [217.6, 302.8] K, max(qˡ): 7.17e-04, Tₛ: 284.9 K, OLR: 293.6 W/m²
[ Info: Iter: 52000, t: 2.242 days, Δt:   4.4s, wall: 18.820 seconds, max|w|:  8.71 m/s, T: [217.6, 302.5] K, max(qˡ): 7.70e-04, Tₛ: 289.2 K, OLR: 297.4 W/m²
[ Info: Iter: 53000, t: 2.290 days, Δt:   4.3s, wall: 19.238 seconds, max|w|: 10.32 m/s, T: [217.7, 302.9] K, max(qˡ): 1.01e-03, Tₛ: 294.6 K, OLR: 303.1 W/m²
[ Info: Iter: 54000, t: 2.333 days, Δt:   3.0s, wall: 18.405 seconds, max|w|: 13.75 m/s, T: [217.5, 302.8] K, max(qˡ): 2.52e-03, Tₛ: 299.9 K, OLR: 307.9 W/m²
[ Info: Iter: 55000, t: 2.372 days, Δt:   3.5s, wall: 18.355 seconds, max|w|:  8.18 m/s, T: [217.6, 303.8] K, max(qˡ): 1.17e-03, Tₛ: 304.8 K, OLR: 313.8 W/m²
[ Info: Iter: 56000, t: 2.411 days, Δt:   3.6s, wall: 18.017 seconds, max|w|: 10.57 m/s, T: [217.7, 302.5] K, max(qˡ): 5.34e-04, Tₛ: 309.4 K, OLR: 319.7 W/m²
[ Info: Iter: 57000, t: 2.451 days, Δt:   4.1s, wall: 18.037 seconds, max|w|: 11.62 m/s, T: [217.3, 304.4] K, max(qˡ): 1.96e-03, Tₛ: 313.5 K, OLR: 323.3 W/m²
[ Info: Iter: 58000, t: 2.492 days, Δt:   2.7s, wall: 18.393 seconds, max|w|: 17.23 m/s, T: [217.1, 302.3] K, max(qˡ): 2.87e-03, Tₛ: 316.8 K, OLR: 325.0 W/m²
[ Info: Iter: 59000, t: 2.525 days, Δt:   3.3s, wall: 17.762 seconds, max|w|: 14.40 m/s, T: [217.6, 304.6] K, max(qˡ): 3.81e-03, Tₛ: 318.7 K, OLR: 328.4 W/m²
[ Info: Iter: 60000, t: 2.561 days, Δt:   2.7s, wall: 17.656 seconds, max|w|: 16.25 m/s, T: [217.6, 306.2] K, max(qˡ): 8.24e-06, Tₛ: 319.8 K, OLR: 330.4 W/m²
[ Info: Iter: 61000, t: 2.597 days, Δt:   2.9s, wall: 17.642 seconds, max|w|: 10.06 m/s, T: [217.7, 304.5] K, max(qˡ): 4.98e-04, Tₛ: 319.9 K, OLR: 328.7 W/m²
[ Info: Iter: 62000, t: 2.632 days, Δt:   2.3s, wall: 18.066 seconds, max|w|: 15.26 m/s, T: [217.7, 304.2] K, max(qˡ): 4.85e-03, Tₛ: 319.1 K, OLR: 327.0 W/m²
[ Info: Iter: 63000, t: 2.663 days, Δt:   2.9s, wall: 17.002 seconds, max|w|:  6.63 m/s, T: [217.6, 304.4] K, max(qˡ): 0.00e+00, Tₛ: 317.5 K, OLR: 325.4 W/m²
[ Info: Iter: 64000, t: 2.704 days, Δt:   4.4s, wall: 18.473 seconds, max|w|:  9.29 m/s, T: [217.5, 303.6] K, max(qˡ): 1.65e-03, Tₛ: 314.5 K, OLR: 322.0 W/m²
[ Info: Iter: 65000, t: 2.749 days, Δt:   3.5s, wall: 18.864 seconds, max|w|: 11.62 m/s, T: [217.6, 305.0] K, max(qˡ): 2.51e-03, Tₛ: 310.1 K, OLR: 315.7 W/m²
[ Info: Iter: 66000, t: 2.798 days, Δt:   4.3s, wall: 19.227 seconds, max|w|:  9.46 m/s, T: [217.5, 303.1] K, max(qˡ): 4.30e-04, Tₛ: 304.4 K, OLR: 309.4 W/m²
[ Info: Iter: 67000, t: 2.842 days, Δt:   4.2s, wall: 18.819 seconds, max|w|: 11.53 m/s, T: [217.7, 303.4] K, max(qˡ): 1.98e-03, Tₛ: 298.9 K, OLR: 302.2 W/m²
[ Info: Iter: 68000, t: 2.882 days, Δt:   2.7s, wall: 18.054 seconds, max|w|: 14.94 m/s, T: [217.5, 304.5] K, max(qˡ): 3.73e-03, Tₛ: 294.0 K, OLR: 298.3 W/m²
[ Info: Iter: 69000, t: 2.923 days, Δt:   4.2s, wall: 18.384 seconds, max|w|:  9.80 m/s, T: [217.1, 301.7] K, max(qˡ): 1.72e-03, Tₛ: 289.3 K, OLR: 293.1 W/m²
[ Info: Iter: 70000, t: 2.967 days, Δt:   3.8s, wall: 13.021 seconds, max|w|: 13.39 m/s, T: [217.3, 302.1] K, max(qˡ): 2.97e-03, Tₛ: 285.1 K, OLR: 289.4 W/m²
[ Info: Simulation is stopping after running for 22.351 minutes.
[ Info: Simulation time 3 days equals or exceeds stop time 3 days.
[ Info: Simulation completed!

Hovmöller diagrams of mean profile evolution

Time–height Hovmöller diagrams show how the horizontally-averaged temperature, cloud liquid water, and radiative flux divergence evolve over two diurnal cycles. The vertical axis is zoomed to the lowest 6 km where the convective dynamics live.

Tts  = FieldTimeSeries(averages_filename, "T")
qᵛts = FieldTimeSeries(averages_filename, "qᵛ")
Fᴿts = FieldTimeSeries(averages_filename, "Fᴿ")

times = Tts.times
Nt = length(times)

z = Oceananigans.Grids.znodes(Tts.grid, Center())
t_hours = times ./ 3600
0.0:1.0:72.0

Build time–height matrices from the averaged profiles.

T_zt  = permutedims(interior(Tts,  1, 1, :, :))
T_zt .-= mean(T_zt, dims=1)
qᵛ_zt = permutedims(interior(qᵛts, 1, 1, :, :))
Fᴿ_zt = permutedims(interior(Fᴿts, 1, 1, :, :))

fig = Figure(size=(900, 800), fontsize=14)

zmax = 12000

axT  = Axis(fig[1, 1]; ylabel="z (m)", title="Temperature anomaly (K)",
            limits=((t_hours[1], t_hours[end]), (0, zmax)))
axqᵛ = Axis(fig[2, 1]; ylabel="z (m)", title="Specific humidity (kg/kg)",
            limits=((t_hours[1], t_hours[end]), (0, zmax)))
axF  = Axis(fig[3, 1]; ylabel="z (m)", xlabel="time (hours)", title="Radiative flux divergence (W/m³)",
            limits=((t_hours[1], t_hours[end]), (0, zmax)))

hidexdecorations!(axT;  grid=false)
hidexdecorations!(axqᵛ; grid=false)

hmT  = heatmap!(axT,  t_hours, z, T_zt;  colormap=:balance, colorrange=(-2, 2))
hmqᵛ = heatmap!(axqᵛ, t_hours, z, qᵛ_zt; colormap=:dense, colorrange=(0, 1e-2))
hmF  = heatmap!(axF,  t_hours, z, Fᴿ_zt; colormap=:balance, colorrange=(-0.1, 0.1))

Colorbar(fig[1, 2], hmT;  label="T′ (K)")
Colorbar(fig[2, 2], hmqᵛ; label="qᵛ (kg/kg)")
Colorbar(fig[3, 2], hmF;  label="Fᴿ (W/m³)")

fig

Animation of cloud structure

We animate xz slices of vertical velocity and cloud liquid water, zoomed to the lowest 5 km where the convective dynamics and clouds live.

wts  = FieldTimeSeries(slices_filename, "w")
qᵛts = FieldTimeSeries(slices_filename, "qᵛ")

times = wts.times
Nt = length(times)

wlim  = maximum(abs, wts) / 4
qᵛlim = maximum(qᵛts) / 2

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

n = Observable(Nt)
axw  = Axis(fig[1, 1]; xlabel="x (km)", ylabel="z (km)", title="w (m/s)", limits=(nothing, (0, 5e3)))
axqᵛ = Axis(fig[1, 2]; xlabel="x (km)", ylabel="z (km)", title="qᵛ (kg/kg)", limits=(nothing, (0, 5e3)))

title = @lift "Diurnal Radiative Convection at t = " * prettytime(times[$n])
fig[0, :] = Label(fig, title, fontsize=16, tellwidth=false)

wn  = @lift wts[$n]
qᵛn = @lift qᵛts[$n]

hmw  = heatmap!(axw,  wn;  colormap=:balance, colorrange=(-wlim, wlim))
hmqᵛ = heatmap!(axqᵛ, qᵛn; colormap=:dense,   colorrange=(0, qᵛlim))

Colorbar(fig[2, 1], hmw;  vertical=false, label="w (m/s)")
Colorbar(fig[2, 2], hmqᵛ; vertical=false, label="qᵛ (g/kg)")

hideydecorations!(axqᵛ; grid=false)

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


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.7
  [13f3f980] CairoMakie v0.15.8
  [6a9e3e04] CloudMicrophysics v0.31.12
  [e30172f5] Documenter v1.17.0
  [daee34ce] DocumenterCitations v1.4.1
  [98b081ad] Literate v2.21.0
  [85f8d34a] NCDatasets v0.14.12
  [9e8cae18] Oceananigans v0.105.2
  [a01a1ee8] RRTMGP v0.21.7
  [b77e0a4c] InteractiveUtils v1.11.0
  [44cfe95a] Pkg v1.12.1
  [9a3f8284] Random v1.11.0

This page was generated using Literate.jl.