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)
if CUDA.functional()
CUDA.seed!(2025)
endGrid
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.0Reference 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 materializedBackground 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 minimum128×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.0Radiation 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
├── 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{south_and_north::Nothing, west_and_east::Nothing, 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!)}}, @NamedTuple{south_and_north::Tuple{Nothing, Nothing}, west_and_east::Tuple{Nothing, Nothing}, bottom_and_top::Tuple{Oceananigans.BoundaryConditions.NoFluxBoundaryCondition, Oceananigans.BoundaryConditions.NoFluxBoundaryCondition}}}, 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{south_and_north::Nothing, west_and_east::Nothing, 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!)}}, @NamedTuple{south_and_north::Tuple{Nothing, Nothing}, west_and_east::Tuple{Nothing, Nothing}, bottom_and_top::Tuple{Oceananigans.BoundaryConditions.BoundaryCondition{Oceananigans.BoundaryConditions.Value, Float32}, Oceananigans.BoundaryConditions.NoFluxBoundaryCondition}}}, 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ᵉ: Centered(order=2)
├── forcing: ρe=>DiscreteForcing
├── tracers: ()
├── coriolis: FPlane{Float32}(f=3.77468e-5)
└── microphysics: SaturationAdjustmentInitial 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ᵗ = specific_prognostic_moisture(model)
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.96983 – 300.16602 K
[ Info: Initial qᵗ range: 0.0 – 12.176381 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)
Oceananigans.Diagnostics.erroring_NaNChecker!(simulation)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: 2.282 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 (24.426 seconds)
[ Info: Executing initial time step...
[ Info: ... initial time step complete (5.705 seconds).
[ Info: Iter: 1000, t: 1.455 hours, Δt: 5.8s, wall: 24.888 seconds, max|w|: 7.69 m/s, T: [217.9, 302.3] K, max(qˡ): 1.96e-03, Tₛ: 280.2 K, OLR: 331.8 W/m²
[ Info: Iter: 2000, t: 2.842 hours, Δt: 3.9s, wall: 15.904 seconds, max|w|: 24.44 m/s, T: [211.9, 302.9] K, max(qˡ): 8.38e-03, Tₛ: 280.5 K, OLR: 183.7 W/m²
[ Info: Iter: 3000, t: 4.114 hours, Δt: 4.3s, wall: 20.225 seconds, max|w|: 15.88 m/s, T: [216.3, 303.3] K, max(qˡ): 7.32e-03, Tₛ: 283.0 K, OLR: 122.4 W/m²
[ Info: Iter: 4000, t: 5.383 hours, Δt: 5.4s, wall: 20.153 seconds, max|w|: 6.18 m/s, T: [216.9, 304.1] K, max(qˡ): 7.87e-03, Tₛ: 287.3 K, OLR: 121.9 W/m²
[ Info: Iter: 5000, t: 6.643 hours, Δt: 4.3s, wall: 20.063 seconds, max|w|: 10.40 m/s, T: [217.2, 306.2] K, max(qˡ): 8.07e-03, Tₛ: 293.0 K, OLR: 122.4 W/m²
[ Info: Iter: 6000, t: 7.772 hours, Δt: 3.3s, wall: 19.549 seconds, max|w|: 12.74 m/s, T: [218.0, 305.0] K, max(qˡ): 7.89e-03, Tₛ: 298.8 K, OLR: 124.3 W/m²
[ Info: Iter: 7000, t: 8.644 hours, Δt: 3.8s, wall: 17.934 seconds, max|w|: 12.04 m/s, T: [218.3, 306.3] K, max(qˡ): 8.24e-03, Tₛ: 303.4 K, OLR: 125.4 W/m²
[ Info: Iter: 8000, t: 9.613 hours, Δt: 3.2s, wall: 18.691 seconds, max|w|: 13.39 m/s, T: [219.1, 305.7] K, max(qˡ): 8.12e-03, Tₛ: 308.2 K, OLR: 126.5 W/m²
[ Info: Iter: 9000, t: 10.434 hours, Δt: 3.8s, wall: 17.885 seconds, max|w|: 9.75 m/s, T: [219.3, 306.7] K, max(qˡ): 8.20e-03, Tₛ: 311.9 K, OLR: 127.1 W/m²
[ Info: Iter: 10000, t: 11.449 hours, Δt: 3.2s, wall: 18.052 seconds, max|w|: 12.34 m/s, T: [219.6, 308.4] K, max(qˡ): 8.45e-03, Tₛ: 315.7 K, OLR: 127.5 W/m²
[ Info: Iter: 11000, t: 12.400 hours, Δt: 4.0s, wall: 18.339 seconds, max|w|: 9.55 m/s, T: [219.1, 306.5] K, max(qˡ): 8.45e-03, Tₛ: 318.3 K, OLR: 127.7 W/m²
[ Info: Iter: 12000, t: 13.401 hours, Δt: 4.2s, wall: 18.672 seconds, max|w|: 9.89 m/s, T: [219.4, 306.5] K, max(qˡ): 8.58e-03, Tₛ: 319.8 K, OLR: 127.5 W/m²
[ Info: Iter: 13000, t: 14.438 hours, Δt: 3.4s, wall: 19.091 seconds, max|w|: 11.97 m/s, T: [219.1, 305.4] K, max(qˡ): 8.65e-03, Tₛ: 319.9 K, OLR: 127.0 W/m²
[ Info: Iter: 14000, t: 15.418 hours, Δt: 4.1s, wall: 18.675 seconds, max|w|: 10.14 m/s, T: [218.7, 305.8] K, max(qˡ): 9.13e-03, Tₛ: 318.6 K, OLR: 126.2 W/m²
[ Info: Iter: 15000, t: 16.386 hours, Δt: 3.8s, wall: 18.296 seconds, max|w|: 10.17 m/s, T: [218.0, 305.4] K, max(qˡ): 8.46e-03, Tₛ: 316.2 K, OLR: 125.1 W/m²
[ Info: Iter: 16000, t: 17.257 hours, Δt: 2.5s, wall: 18.285 seconds, max|w|: 17.62 m/s, T: [217.6, 307.8] K, max(qˡ): 8.60e-03, Tₛ: 313.2 K, OLR: 123.8 W/m²
[ Info: Iter: 17000, t: 18.208 hours, Δt: 4.4s, wall: 18.011 seconds, max|w|: 10.92 m/s, T: [217.2, 306.2] K, max(qˡ): 8.51e-03, Tₛ: 309.0 K, OLR: 122.3 W/m²
[ Info: Iter: 18000, t: 19.181 hours, Δt: 4.3s, wall: 19.228 seconds, max|w|: 10.46 m/s, T: [216.7, 306.1] K, max(qˡ): 7.54e-03, Tₛ: 304.3 K, OLR: 122.0 W/m²
[ Info: Iter: 19000, t: 20.290 hours, Δt: 5.2s, wall: 19.033 seconds, max|w|: 6.35 m/s, T: [216.3, 304.9] K, max(qˡ): 8.59e-03, Tₛ: 298.5 K, OLR: 121.9 W/m²
[ Info: Iter: 20000, t: 21.587 hours, Δt: 4.7s, wall: 20.269 seconds, max|w|: 8.92 m/s, T: [217.0, 304.1] K, max(qˡ): 8.14e-03, Tₛ: 291.9 K, OLR: 121.9 W/m²
[ Info: Iter: 21000, t: 22.706 hours, Δt: 3.6s, wall: 19.099 seconds, max|w|: 8.27 m/s, T: [216.7, 305.0] K, max(qˡ): 7.77e-03, Tₛ: 287.0 K, OLR: 121.9 W/m²
[ Info: Iter: 22000, t: 23.950 hours, Δt: 3.9s, wall: 19.464 seconds, max|w|: 6.89 m/s, T: [216.9, 304.1] K, max(qˡ): 7.50e-03, Tₛ: 282.8 K, OLR: 121.9 W/m²
[ Info: Iter: 23000, t: 1.044 days, Δt: 3.8s, wall: 18.925 seconds, max|w|: 11.63 m/s, T: [217.2, 304.8] K, max(qˡ): 7.86e-03, Tₛ: 280.6 K, OLR: 121.9 W/m²
[ Info: Iter: 24000, t: 1.080 days, Δt: 3.9s, wall: 17.873 seconds, max|w|: 8.44 m/s, T: [216.8, 306.2] K, max(qˡ): 8.27e-03, Tₛ: 280.0 K, OLR: 122.0 W/m²
[ Info: Iter: 25000, t: 1.121 days, Δt: 4.5s, wall: 18.616 seconds, max|w|: 9.57 m/s, T: [216.3, 304.5] K, max(qˡ): 7.90e-03, Tₛ: 280.6 K, OLR: 121.9 W/m²
[ Info: Iter: 26000, t: 1.159 days, Δt: 2.9s, wall: 18.242 seconds, max|w|: 15.64 m/s, T: [217.1, 306.9] K, max(qˡ): 7.36e-03, Tₛ: 282.2 K, OLR: 122.0 W/m²
[ Info: Iter: 27000, t: 1.199 days, Δt: 4.5s, wall: 18.638 seconds, max|w|: 10.50 m/s, T: [216.9, 306.6] K, max(qˡ): 7.42e-03, Tₛ: 285.0 K, OLR: 121.8 W/m²
[ Info: Iter: 28000, t: 1.242 days, Δt: 4.0s, wall: 18.701 seconds, max|w|: 11.15 m/s, T: [216.8, 303.7] K, max(qˡ): 7.46e-03, Tₛ: 289.1 K, OLR: 121.9 W/m²
[ Info: Iter: 29000, t: 1.283 days, Δt: 4.0s, wall: 18.707 seconds, max|w|: 11.56 m/s, T: [217.3, 305.5] K, max(qˡ): 7.31e-03, Tₛ: 293.8 K, OLR: 122.6 W/m²
[ Info: Iter: 30000, t: 1.325 days, Δt: 4.2s, wall: 18.682 seconds, max|w|: 10.04 m/s, T: [218.2, 304.5] K, max(qˡ): 7.92e-03, Tₛ: 298.9 K, OLR: 124.2 W/m²
[ Info: Iter: 31000, t: 1.363 days, Δt: 3.5s, wall: 18.309 seconds, max|w|: 11.09 m/s, T: [218.8, 307.8] K, max(qˡ): 7.23e-03, Tₛ: 303.7 K, OLR: 125.5 W/m²
[ Info: Iter: 32000, t: 1.404 days, Δt: 2.9s, wall: 18.651 seconds, max|w|: 13.53 m/s, T: [219.2, 303.8] K, max(qˡ): 7.37e-03, Tₛ: 308.6 K, OLR: 126.5 W/m²
[ Info: Iter: 33000, t: 1.446 days, Δt: 3.1s, wall: 18.670 seconds, max|w|: 12.89 m/s, T: [219.2, 306.0] K, max(qˡ): 7.39e-03, Tₛ: 313.0 K, OLR: 127.2 W/m²
[ Info: Iter: 34000, t: 1.492 days, Δt: 3.6s, wall: 19.890 seconds, max|w|: 10.66 m/s, T: [219.6, 308.7] K, max(qˡ): 7.02e-03, Tₛ: 316.8 K, OLR: 127.6 W/m²
[ Info: Iter: 35000, t: 1.533 days, Δt: 3.7s, wall: 18.078 seconds, max|w|: 12.60 m/s, T: [219.6, 305.6] K, max(qˡ): 7.23e-03, Tₛ: 319.0 K, OLR: 127.6 W/m²
[ Info: Iter: 36000, t: 1.576 days, Δt: 3.6s, wall: 18.657 seconds, max|w|: 15.04 m/s, T: [219.5, 305.1] K, max(qˡ): 6.94e-03, Tₛ: 320.0 K, OLR: 127.3 W/m²
[ Info: Iter: 37000, t: 1.615 days, Δt: 3.5s, wall: 18.557 seconds, max|w|: 11.34 m/s, T: [219.1, 306.5] K, max(qˡ): 7.02e-03, Tₛ: 319.6 K, OLR: 126.8 W/m²
[ Info: Iter: 38000, t: 1.656 days, Δt: 4.6s, wall: 18.245 seconds, max|w|: 9.85 m/s, T: [218.6, 305.3] K, max(qˡ): 6.82e-03, Tₛ: 317.9 K, OLR: 125.9 W/m²
[ Info: Iter: 39000, t: 1.691 days, Δt: 2.1s, wall: 18.250 seconds, max|w|: 20.84 m/s, T: [218.3, 304.5] K, max(qˡ): 7.01e-03, Tₛ: 315.6 K, OLR: 124.8 W/m²
[ Info: Iter: 40000, t: 1.727 days, Δt: 3.7s, wall: 17.922 seconds, max|w|: 9.22 m/s, T: [217.9, 304.2] K, max(qˡ): 7.12e-03, Tₛ: 312.4 K, OLR: 123.6 W/m²
[ Info: Iter: 41000, t: 1.768 days, Δt: 3.9s, wall: 18.738 seconds, max|w|: 7.92 m/s, T: [217.3, 305.0] K, max(qˡ): 7.12e-03, Tₛ: 308.0 K, OLR: 122.0 W/m²
[ Info: Iter: 42000, t: 1.805 days, Δt: 4.2s, wall: 17.913 seconds, max|w|: 8.18 m/s, T: [217.2, 306.3] K, max(qˡ): 7.00e-03, Tₛ: 303.5 K, OLR: 121.9 W/m²
[ Info: Iter: 43000, t: 1.851 days, Δt: 3.5s, wall: 19.451 seconds, max|w|: 9.07 m/s, T: [217.3, 304.4] K, max(qˡ): 7.36e-03, Tₛ: 297.7 K, OLR: 121.9 W/m²
[ Info: Iter: 44000, t: 1.890 days, Δt: 3.9s, wall: 18.335 seconds, max|w|: 10.68 m/s, T: [217.1, 307.2] K, max(qˡ): 7.28e-03, Tₛ: 293.0 K, OLR: 121.9 W/m²
[ Info: Iter: 45000, t: 1.931 days, Δt: 2.4s, wall: 18.757 seconds, max|w|: 16.41 m/s, T: [217.3, 304.2] K, max(qˡ): 7.27e-03, Tₛ: 288.5 K, OLR: 122.0 W/m²
[ Info: Iter: 46000, t: 1.969 days, Δt: 3.6s, wall: 17.947 seconds, max|w|: 7.55 m/s, T: [216.8, 304.9] K, max(qˡ): 6.63e-03, Tₛ: 285.0 K, OLR: 121.9 W/m²
[ Info: Iter: 47000, t: 2.017 days, Δt: 4.5s, wall: 18.908 seconds, max|w|: 9.89 m/s, T: [217.3, 307.4] K, max(qˡ): 6.93e-03, Tₛ: 281.7 K, OLR: 121.9 W/m²
[ Info: Iter: 48000, t: 2.062 days, Δt: 4.0s, wall: 19.064 seconds, max|w|: 10.12 m/s, T: [217.2, 303.7] K, max(qˡ): 6.86e-03, Tₛ: 280.2 K, OLR: 121.9 W/m²
[ Info: Iter: 49000, t: 2.106 days, Δt: 3.6s, wall: 18.538 seconds, max|w|: 12.00 m/s, T: [217.1, 304.2] K, max(qˡ): 6.90e-03, Tₛ: 280.2 K, OLR: 121.9 W/m²
[ Info: Iter: 50000, t: 2.147 days, Δt: 4.4s, wall: 19.503 seconds, max|w|: 12.47 m/s, T: [217.2, 304.4] K, max(qˡ): 7.01e-03, Tₛ: 281.6 K, OLR: 121.9 W/m²
[ Info: Iter: 51000, t: 2.186 days, Δt: 3.6s, wall: 18.243 seconds, max|w|: 10.86 m/s, T: [217.2, 304.8] K, max(qˡ): 6.96e-03, Tₛ: 284.1 K, OLR: 121.9 W/m²
[ Info: Iter: 52000, t: 2.225 days, Δt: 3.1s, wall: 18.285 seconds, max|w|: 17.34 m/s, T: [217.0, 306.8] K, max(qˡ): 6.64e-03, Tₛ: 287.5 K, OLR: 121.9 W/m²
[ Info: Iter: 53000, t: 2.265 days, Δt: 4.4s, wall: 18.707 seconds, max|w|: 12.45 m/s, T: [217.2, 304.1] K, max(qˡ): 6.91e-03, Tₛ: 291.6 K, OLR: 122.0 W/m²
[ Info: Iter: 54000, t: 2.309 days, Δt: 3.8s, wall: 19.107 seconds, max|w|: 11.13 m/s, T: [217.9, 305.4] K, max(qˡ): 6.54e-03, Tₛ: 297.0 K, OLR: 123.7 W/m²
[ Info: Iter: 55000, t: 2.353 days, Δt: 3.8s, wall: 18.717 seconds, max|w|: 10.24 m/s, T: [218.7, 306.2] K, max(qˡ): 6.40e-03, Tₛ: 302.4 K, OLR: 125.2 W/m²
[ Info: Iter: 56000, t: 2.398 days, Δt: 3.1s, wall: 19.072 seconds, max|w|: 14.31 m/s, T: [219.1, 306.7] K, max(qˡ): 6.44e-03, Tₛ: 307.9 K, OLR: 126.4 W/m²
[ Info: Iter: 57000, t: 2.443 days, Δt: 4.6s, wall: 19.080 seconds, max|w|: 8.83 m/s, T: [219.4, 307.0] K, max(qˡ): 6.46e-03, Tₛ: 312.7 K, OLR: 127.2 W/m²
[ Info: Iter: 58000, t: 2.486 days, Δt: 4.2s, wall: 18.635 seconds, max|w|: 8.33 m/s, T: [219.8, 305.5] K, max(qˡ): 6.42e-03, Tₛ: 316.3 K, OLR: 127.6 W/m²
[ Info: Iter: 59000, t: 2.527 days, Δt: 3.2s, wall: 18.136 seconds, max|w|: 14.69 m/s, T: [219.6, 306.0] K, max(qˡ): 6.42e-03, Tₛ: 318.8 K, OLR: 127.6 W/m²
[ Info: Iter: 60000, t: 2.570 days, Δt: 4.0s, wall: 19.060 seconds, max|w|: 11.13 m/s, T: [219.5, 305.9] K, max(qˡ): 6.51e-03, Tₛ: 319.9 K, OLR: 127.4 W/m²
[ Info: Iter: 61000, t: 2.614 days, Δt: 3.7s, wall: 18.665 seconds, max|w|: 13.09 m/s, T: [219.3, 307.3] K, max(qˡ): 6.46e-03, Tₛ: 319.6 K, OLR: 126.8 W/m²
[ Info: Iter: 62000, t: 2.658 days, Δt: 5.0s, wall: 19.044 seconds, max|w|: 9.91 m/s, T: [219.0, 306.4] K, max(qˡ): 6.62e-03, Tₛ: 317.9 K, OLR: 125.8 W/m²
[ Info: Iter: 63000, t: 2.700 days, Δt: 3.6s, wall: 18.676 seconds, max|w|: 9.58 m/s, T: [218.2, 306.8] K, max(qˡ): 6.79e-03, Tₛ: 314.8 K, OLR: 124.6 W/m²
[ Info: Iter: 64000, t: 2.743 days, Δt: 3.7s, wall: 18.720 seconds, max|w|: 9.58 m/s, T: [217.6, 306.2] K, max(qˡ): 6.61e-03, Tₛ: 310.8 K, OLR: 123.0 W/m²
[ Info: Iter: 65000, t: 2.782 days, Δt: 3.3s, wall: 18.661 seconds, max|w|: 16.22 m/s, T: [217.2, 308.9] K, max(qˡ): 6.37e-03, Tₛ: 306.3 K, OLR: 121.9 W/m²
[ Info: Iter: 66000, t: 2.814 days, Δt: 2.7s, wall: 18.288 seconds, max|w|: 7.03 m/s, T: [217.3, 305.7] K, max(qˡ): 6.58e-03, Tₛ: 302.4 K, OLR: 122.0 W/m²
[ Info: Iter: 67000, t: 2.856 days, Δt: 3.5s, wall: 18.681 seconds, max|w|: 11.39 m/s, T: [217.0, 305.1] K, max(qˡ): 6.16e-03, Tₛ: 297.2 K, OLR: 121.8 W/m²
[ Info: Iter: 68000, t: 2.899 days, Δt: 4.4s, wall: 18.689 seconds, max|w|: 12.01 m/s, T: [217.2, 305.5] K, max(qˡ): 6.42e-03, Tₛ: 292.0 K, OLR: 121.9 W/m²
[ Info: Iter: 69000, t: 2.939 days, Δt: 2.8s, wall: 18.742 seconds, max|w|: 16.43 m/s, T: [217.2, 304.6] K, max(qˡ): 6.31e-03, Tₛ: 287.7 K, OLR: 121.9 W/m²
[ Info: Iter: 70000, t: 2.981 days, Δt: 4.3s, wall: 18.689 seconds, max|w|: 12.78 m/s, T: [217.3, 304.6] K, max(qˡ): 6.36e-03, Tₛ: 284.0 K, OLR: 121.9 W/m²
[ Info: Simulation is stopping after running for 22.364 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 ./ 36000.0:1.0:72.0Build 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³)")
figAnimation 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
endJulia version and environment information
This example was executed with the following version of Julia:
using InteractiveUtils: versioninfo
versioninfo()Julia Version 1.12.5
Commit 5fe89b8ddc1 (2026-02-09 16:05 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_PKG_SERVER_REGISTRY_PREFERENCE = eager
JULIA_VERSION = 1.12.5
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.4.6 `.`
[052768ef] CUDA v5.11.0
[13f3f980] CairoMakie v0.15.9
⌅ [6a9e3e04] CloudMicrophysics v0.32.0
[e30172f5] Documenter v1.17.0
[daee34ce] DocumenterCitations v1.4.1
[98b081ad] Literate v2.21.0
[85f8d34a] NCDatasets v0.14.15
[9e8cae18] Oceananigans v0.106.4
[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.