Prescribed sea surface temperature convection
This example simulates moist convection driven by a prescribed sea surface temperature (SST). The simulation models the atmospheric response to a horizontally-varying SST pattern, a fundamental problem in atmosphere-ocean interaction studies. The setup is representative of convection over oceanic fronts or sea surface temperature gradients, where differential heating drives organized atmospheric circulations.
The simulation uses bulk aerodynamic formulas to compute surface fluxes of momentum, sensible heat, and latent heat based on bulk transfer coefficients. This approach parameterizes the complex turbulent exchange processes in the surface layer using simple drag law formulations that relate fluxes to the difference between surface and near-surface atmospheric properties.
The model uses warm-phase saturation adjustment microphysics with liquid-ice potential temperature thermodynamics. Saturation adjustment instantly condenses or evaporates water vapor to maintain thermodynamic equilibrium, providing a simple yet effective representation of cloud processes in moist convection.
using Breeze
using Oceananigans
using Oceananigans.Units
using Oceananigans.Models: BoundaryConditionOperation
using Printf
using CairoMakieGrid setup
We use a 2D domain (x-z plane) with periodic horizontal boundaries and a bounded vertical domain. The horizontal periodicity allows convective cells to develop and interact without artificial boundary effects. The domain extends 20 km horizontally and 10 km vertically.
The grid resolution of 128 points in each direction provides approximately 156 m horizontal and 78 m vertical resolution, sufficient to resolve the energy-containing scales of convective turbulence while remaining computationally tractable for this demonstration.
grid = RectilinearGrid(size = (128, 128), halo = (5, 5),
x = (-10kilometers, 10kilometers),
z = (0, 10kilometers),
topology = (Periodic, Flat, Bounded))128×1×128 RectilinearGrid{Float64, Periodic, Flat, Bounded} on CPU with 5×0×5 halo
├── Periodic x ∈ [-10000.0, 10000.0) regularly spaced with Δx=156.25
├── Flat y
└── Bounded z ∈ [0.0, 10000.0] regularly spaced with Δz=78.125Model formulation
We create an AtmosphereModel with warm-phase saturation adjustment microphysics and liquid-ice potential temperature thermodynamics. The anelastic formulation filters acoustic waves while retaining the essential dynamics of deep convection, allowing larger time steps than a fully compressible model.
The reference state defines the background thermodynamic profile against which perturbations evolve. We use a base pressure p₀ = 101325 Pa (standard sea level pressure) and reference potential temperature θ₀ = 285 K, representing a relatively cool maritime atmosphere.
p₀, θ₀ = 101325, 285 # Pa, K
constants = ThermodynamicConstants()
reference_state = ReferenceState(grid, constants; surface_pressure=p₀, potential_temperature=θ₀)
dynamics = AnelasticDynamics(reference_state)AnelasticDynamics(p₀=101325.0, θ₀=285.0)
└── pressure_anomaly: not materializedThe microphysics scheme uses saturation adjustment to maintain thermodynamic equilibrium. The WarmPhaseEquilibrium option considers only liquid water and vapor, appropriate for warm convection where ice processes are negligible.
microphysics = SaturationAdjustment(equilibrium = WarmPhaseEquilibrium())Breeze.Microphysics.SaturationAdjustment{Breeze.Thermodynamics.WarmPhaseEquilibrium, Float64}(0.001, Inf, Breeze.Thermodynamics.WarmPhaseEquilibrium())We use high-order WENO advection schemes to accurately represent the sharp gradients that develop in convective flows. WENO (Weighted Essentially Non-Oscillatory; Shu (2009)) schemes provide excellent shock-capturing properties while maintaining high accuracy in smooth regions.
momentum_advection = WENO(order=9)
scalar_advection = WENO(order=5)WENO{3, Float64, Float32}(order=5)
├── buffer_scheme: WENO{2, Float64, Float32}(order=3)
│ └── buffer_scheme: Centered(order=2)
└── advecting_velocity_scheme: Centered(order=4)Boundary conditions
Breeze provides abstractions for specifying bulk surface fluxes. The BulkDrag, BulkSensibleHeatFlux, and BulkVaporFlux boundary conditions compute fluxes of momentum, temperature density (proportional to sensible heat flux), and moisture density according to bulk aerodynamic formulae that relate turbulent fluxes to the difference between atmosphere properties, surface properties, and the differential motion of the air and surface,
\[τˣ = - Cᴰ |U| ρu, \quad Jᵀ = - ρ₀ Cᵀ |U| (θ - θ₀), \quad Jᵛ = - ρ₀ Cᵛ |U| (qᵗ - qᵛ₀),\]
where $|U|$ is "total" the differential wind speed (including gustiness), $Cᴰ, Cᵀ, Cᵛ$ are transfer coefficients, and $θ₀, qᵛ₀$ are the surface temperature and surface specific humidity, which for wet surfaces is presumed to be the saturation specific humidity over a planar liquid surface computed at the surface temperature. $τˣ$ is the surface momentum flux, $Jᵀ$ is the surface temperature density flux, and $Jᵛ$ is the surface moisture density flux. The surface density density $ρ₀$ is computed from the model's reference state.
The temperature density flux is proportional to the sensible heat flux,
\[𝒬ᵀ = - ρ₀ cᵖᵐ Cᵀ |U| (θ - θ₀) .\]
where $cᵖᵐ$ is the mixture heat capacity.
We start by defining the drag coefficient and gustiness parameter,
Cᴰ = 1e-3 # Drag coefficient
Uᵍ = 1e-2 # Minimum wind speed (m/s)
ρu_surface_flux = ρv_surface_flux = BulkDrag(coefficient=Cᴰ, gustiness=Uᵍ)FluxBoundaryCondition: BulkDragFunction(direction=Nothing, coefficient=0.001, gustiness=0.01)Sensible heat flux and vapor fluxes
For BulkVaporFlux, the saturation specific humidity is computed from the surface temperature. Surface temperature can be provided as a Field, a Function, or a Number.
In this example, we specify the sea surface temperature as a top hat function i.e. representing a pair of ocean fronts in a periodic domain, with a difference of 4 degrees K,
ΔT = 4 # K
T₀(x) = θ₀ + ΔT / 2 * sign(cos(2π * x / grid.Lx))T₀ (generic function with 1 method)We complete our specification with the sensible heat transfer coefficient and vapor transfer coefficient,
Cᵀ = 1e-3 # Sensible heat transfer coefficient
Cᵛ = 1e-3 # Vapor transfer coefficient0.001and build the flux parameterizations
ρθ_surface_flux = BulkSensibleHeatFlux(coefficient=Cᵀ, gustiness=Uᵍ, surface_temperature=T₀)
ρqᵗ_surface_flux = BulkVaporFlux(coefficient=Cᵛ, gustiness=Uᵍ, surface_temperature=T₀)FluxBoundaryCondition: BulkVaporFluxFunction(coefficient=0.001, gustiness=0.01)We finally assemble all of the boundary conditions,
ρu_bcs = FieldBoundaryConditions(bottom=ρu_surface_flux)
ρv_bcs = FieldBoundaryConditions(bottom=ρv_surface_flux)
ρθ_bcs = FieldBoundaryConditions(bottom=ρθ_surface_flux)
ρqᵗ_bcs = FieldBoundaryConditions(bottom=ρqᵗ_surface_flux)Oceananigans.FieldBoundaryConditions, with boundary conditions
├── west: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── east: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── south: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── north: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── bottom: FluxBoundaryCondition: BulkVaporFluxFunction(coefficient=0.001, gustiness=0.01)
├── top: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
└── immersed: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)Model construction
We assemble the AtmosphereModel with all the components defined above. The model will solve the anelastic equations with the specified advection schemes, microphysics, and boundary conditions.
model = AtmosphereModel(grid; momentum_advection, scalar_advection, microphysics, dynamics,
boundary_conditions = (ρu=ρu_bcs, ρv=ρv_bcs, ρθ=ρθ_bcs, ρqᵗ=ρqᵗ_bcs))AtmosphereModel{CPU, RectilinearGrid}(time = 0 seconds, iteration = 0)
├── grid: 128×1×128 RectilinearGrid{Float64, Periodic, Flat, Bounded} on CPU with 5×0×5 halo
├── dynamics: AnelasticDynamics(p₀=101325.0, θ₀=285.0)
├── formulation: LiquidIcePotentialTemperatureFormulation
├── timestepper: SSPRungeKutta3
├── advection scheme:
│ ├── momentum: WENO{5, Float64, Float32}(order=9)
│ ├── ρθ: WENO{3, Float64, Float32}(order=5)
│ └── ρqᵗ: WENO{3, Float64, Float32}(order=5)
├── tracers: ()
├── coriolis: Nothing
└── microphysics: SaturationAdjustmentInitial conditions
We initialize the model with a uniform potential temperature equal to the reference value, creating a neutrally stratified atmosphere. A small background wind (1 m/s) in the x-direction provides initial momentum for the bulk flux calculations and helps break symmetry.
set!(model, θ=reference_state.potential_temperature, u=1)Simulation setup
We configure the simulation to run for 4 hours with adaptive time stepping. The CFL condition limits the time step to maintain numerical stability, with a target CFL number of 0.7 providing a good balance between efficiency and accuracy.
simulation = Simulation(model, Δt=10, stop_time=4hours)
conjure_time_step_wizard!(simulation, cfl=0.7)Diagnostic fields
We define several diagnostic quantities for analysis and visualization:
- Temperature T: the actual temperature field
- Potential temperature θ: conserved in dry adiabatic processes
- Liquid water content qˡ: mass fraction of cloud liquid water
- Saturation specific humidity qᵛ⁺: maximum water vapor the air can hold
T = model.temperature
θ = liquid_ice_potential_temperature(model)
qˡ = model.microphysical_fields.qˡ
qᵛ⁺ = Breeze.Microphysics.SaturationSpecificHumidity(model)
ρu, ρv, ρw = model.momentum
u, v, w = model.velocities
qᵗ = model.specific_moisture128×1×128 Field{Oceananigans.Grids.Center, Oceananigans.Grids.Center, Oceananigans.Grids.Center} on Oceananigans.Grids.RectilinearGrid on CPU
├── grid: 128×1×128 RectilinearGrid{Float64, Periodic, Flat, Bounded} on CPU with 5×0×5 halo
├── boundary conditions: FieldBoundaryConditions
│ └── west: Periodic, east: Periodic, south: Nothing, north: Nothing, bottom: ZeroFlux, top: ZeroFlux, immersed: Nothing
└── data: 138×1×138 OffsetArray(::Array{Float64, 3}, -4:133, 1:1, -4:133) with eltype Float64 with indices -4:133×1:1×-4:133
└── max=0.0, min=0.0, mean=0.0Surface flux diagnostics
We use Oceananigans' BoundaryConditionOperation to extract the surface flux values from the boundary conditions. These 1D fields (varying only in x) represent the actual flux values applied at the ocean-atmosphere interface.
The surface fluxes are:
- $τˣ$: momentum flux (stress), in kg m⁻¹ s⁻²
- $𝒬ᵀ$: sensible heat flux = cᵖᵐ Jᵀ, in W m⁻²
- $𝒬ᵛ$: latent heat flux = ℒˡ Jᵛ, in W m⁻²
where Jᵀ is the temperature density flux and Jᵛ is the moisture density flux.
# Surface momentum flux
τˣ = BoundaryConditionOperation(ρu, :bottom, model)
# Sensible heat flux: 𝒬ᵀ = cᵖᵈ Jᵀ (using dry air heat capacity as approximation)
ρθ = liquid_ice_potential_temperature_density(model)
cᵖᵈ = constants.dry_air.heat_capacity
Jᵀ = BoundaryConditionOperation(ρθ, :bottom, model)
𝒬ᵀ = cᵖᵈ * Jᵀ
# Latent heat flux: 𝒬ᵛ = ℒˡ Jᵛ (using reference θ₀ for latent heat)
ρqᵗ = model.moisture_density
ℒˡ = Breeze.Thermodynamics.liquid_latent_heat(θ₀, constants)
Jᵛ = BoundaryConditionOperation(ρqᵗ, :bottom, model)
𝒬ᵛ = ℒˡ * JᵛBinaryOperation at (Center, Center, ⋅)
├── grid: 128×1×128 RectilinearGrid{Float64, Periodic, Flat, Bounded} on CPU with 5×0×5 halo
└── tree:
* at (Center, Center, ⋅)
├── 2.47317765e6
└── KernelFunctionOperation at (Center, Center, ⋅)Progress callback
A callback function prints diagnostic information every few iterations, helping monitor the simulation's progress and detect any numerical issues.
function progress(sim)
qᵗ = sim.model.specific_moisture
u, v, w = sim.model.velocities
umax = maximum(abs, u)
vmax = maximum(abs, v)
wmax = maximum(abs, w)
qᵗmin = minimum(qᵗ)
qᵗmax = maximum(qᵗ)
qˡmax = maximum(qˡ)
θmin = minimum(θ)
θmax = maximum(θ)
msg = @sprintf("Iter: %d, t = %s, max|u|: (%.2e, %.2e, %.2e)",
iteration(sim), prettytime(sim), umax, vmax, wmax)
msg *= @sprintf(", extrema(qᵗ): (%.2e, %.2e), max(qˡ): %.2e, extrema(θ): (%.2e, %.2e)",
qᵗmin, qᵗmax, qˡmax, θmin, θmax)
@info msg
return nothing
end
add_callback!(simulation, progress, IterationInterval(100))Output
We save both the full 2D fields and the 1D surface flux fields. We include both native model variables and others like, e.g., the total speed, $\sqrt{u² + w²}$ and the cross-stream vorticity $∂_z u - ∂_x w$. The JLD2 format provides efficient storage with full Julia type preservation.
output_filename = "prescribed_sea_surface_temperature_convection.jld2"
qᵗ = model.specific_moisture
u, v, w, = model.velocities
s = sqrt(u^2 + w^2) # speed
ξ = ∂z(u) - ∂x(w) # cross-stream vorticity
outputs = (; s, ξ, T, θ, qˡ, qᵛ⁺, qᵗ, τˣ, 𝒬ᵀ, 𝒬ᵛ, Σ𝒬=𝒬ᵀ+𝒬ᵛ)
ow = JLD2Writer(model, outputs;
filename = output_filename,
schedule = TimeInterval(2minutes),
overwrite_existing = true)
simulation.output_writers[:jld2] = owJLD2Writer scheduled on TimeInterval(2 minutes):
├── filepath: prescribed_sea_surface_temperature_convection.jld2
├── 11 outputs: (s, ξ, T, θ, qˡ, qᵛ⁺, qᵗ, τˣ, 𝒬ᵀ, 𝒬ᵛ, Σ𝒬)
├── array_type: Array{Float32}
├── including: [:grid]
├── file_splitting: NoFileSplitting
└── file size: 70.7 KiBRun the simulation
@info "Running prescribed SST convection simulation..."
run!(simulation)[ Info: Running prescribed SST convection simulation...
[ Info: Initializing simulation...
[ Info: Iter: 0, t = 0 seconds, max|u|: (1.00e+00, 0.00e+00, 0.00e+00), extrema(qᵗ): (0.00e+00, 0.00e+00), max(qˡ): 0.00e+00, extrema(θ): (2.85e+02, 2.85e+02)
[ Info: ... simulation initialization complete (30.829 seconds)
[ Info: Executing initial time step...
[ Info: ... initial time step complete (4.700 seconds).
[ Info: Iter: 100, t = 27.297 minutes, max|u|: (1.08e+00, 0.00e+00, 2.27e-02), extrema(qᵗ): (-6.72e-06, 2.13e-04), max(qˡ): 0.00e+00, extrema(θ): (2.85e+02, 2.85e+02)
[ Info: Iter: 200, t = 1.433 hours, max|u|: (2.41e+00, 0.00e+00, 9.50e-01), extrema(qᵗ): (-1.52e-04, 7.19e-04), max(qˡ): 0.00e+00, extrema(θ): (2.85e+02, 2.85e+02)
[ Info: Iter: 300, t = 2.158 hours, max|u|: (2.55e+00, 0.00e+00, 2.47e+00), extrema(qᵗ): (-2.80e-04, 1.01e-03), max(qˡ): 0.00e+00, extrema(θ): (2.85e+02, 2.85e+02)
[ Info: Iter: 400, t = 2.700 hours, max|u|: (2.76e+00, 0.00e+00, 1.85e+00), extrema(qᵗ): (-2.36e-04, 1.09e-03), max(qˡ): 0.00e+00, extrema(θ): (2.85e+02, 2.85e+02)
[ Info: Iter: 500, t = 3.227 hours, max|u|: (2.68e+00, 0.00e+00, 2.26e+00), extrema(qᵗ): (-2.58e-04, 1.25e-03), max(qˡ): 0.00e+00, extrema(θ): (2.85e+02, 2.85e+02)
[ Info: Iter: 600, t = 3.750 hours, max|u|: (3.15e+00, 0.00e+00, 2.04e+00), extrema(qᵗ): (-3.71e-04, 1.27e-03), max(qˡ): 2.59e-06, extrema(θ): (2.85e+02, 2.85e+02)
[ Info: Simulation is stopping after running for 1.540 minutes.
[ Info: Simulation time 4 hours equals or exceeds stop time 4 hours.
Visualization
We create animations showing the evolution of the flow fields. The figure displays velocity components (u, w), thermodynamic fields (θ, T), moisture fields (qᵗ, qˡ), and surface fluxes (momentum and heat).
@assert isfile(output_filename) "Output file $(output_filename) not found."
s_ts = FieldTimeSeries(output_filename, "s")
ξ_ts = FieldTimeSeries(output_filename, "ξ")
θ_ts = FieldTimeSeries(output_filename, "θ")
T_ts = FieldTimeSeries(output_filename, "T")
qᵗ_ts = FieldTimeSeries(output_filename, "qᵗ")
qˡ_ts = FieldTimeSeries(output_filename, "qˡ")
τˣ_ts = FieldTimeSeries(output_filename, "τˣ")
𝒬ᵀ_ts = FieldTimeSeries(output_filename, "𝒬ᵀ")
𝒬ᵛ_ts = FieldTimeSeries(output_filename, "𝒬ᵛ")
Σ𝒬_ts = FieldTimeSeries(output_filename, "Σ𝒬")
times = θ_ts.times
Nt = length(θ_ts)
n = Observable(Nt)
sn = @lift s_ts[$n]
ξn = @lift ξ_ts[$n]
θn = @lift θ_ts[$n]
qᵗn = @lift qᵗ_ts[$n]
Tn = @lift T_ts[$n]
qˡn = @lift qˡ_ts[$n]
τˣn = @lift τˣ_ts[$n]
𝒬ᵀn = @lift 𝒬ᵀ_ts[$n]
𝒬ᵛn = @lift 𝒬ᵛ_ts[$n]
Σ𝒬n = @lift Σ𝒬_ts[$n]Observable(128×1×1 Field{Oceananigans.Grids.Center, Oceananigans.Grids.Center, Nothing} reduced over dims = (3,) on Oceananigans.Grids.RectilinearGrid on CPU
├── grid: 128×1×128 RectilinearGrid{Float64, Periodic, Flat, Bounded} on CPU 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(view(::Array{Float64, 4}, :, :, :, 121), -4:133, 1:1, 1:1) with eltype Float64 with indices -4:133×1:1×1:1
└── max=79.2817, min=0.522659, mean=29.385)
Now we are ready to plot.
fig = Figure(size=(800, 1000), fontsize=13)
title = @lift "t = $(prettytime(times[$n]))"
axs = Axis(fig[1, 1], ylabel="z (m)")
axξ = Axis(fig[1, 2])
axθ = Axis(fig[2, 1], ylabel="z (m)")
axq = Axis(fig[2, 2])
axT = Axis(fig[3, 1], ylabel="z (m)")
axqˡ = Axis(fig[3, 2])Axis with 0 plots:
Surface flux plots at bottom
axτ = Axis(fig[4, 1], xlabel="x (m)", ylabel="τˣ (kg m⁻¹ s⁻²)", title="Surface momentum flux")
ax𝒬 = Axis(fig[4, 2], xlabel="x (m)", ylabel="𝒬 (W m⁻²)", title="Surface heat flux (𝒬ᵀ + 𝒬ᵛ)")
fig[0, :] = Label(fig, title, fontsize=22, tellwidth=false)Label()Compute color limits from the full time series
θ_limits = extrema(θ_ts)
T_limits = extrema(T_ts)
s_limits = (0, maximum(s_ts))
ξ_lim = 0.8 * maximum(abs, ξ_ts)
ξ_limits = (-ξ_lim, +ξ_lim)
qᵗ_max = maximum(qᵗ_ts)
qˡ_max = maximum(qˡ_ts)1.8863338482333347e-5Flux limits
τˣ_max = max(abs(minimum(τˣ_ts)), abs(maximum(τˣ_ts)))
𝒬_min = min(minimum(𝒬ᵀ_ts), minimum(𝒬ᵛ_ts), minimum(Σ𝒬_ts))
𝒬_max = max(maximum(𝒬ᵀ_ts), maximum(𝒬ᵛ_ts), maximum(Σ𝒬_ts))
hms = heatmap!(axs, sn, colorrange=s_limits, colormap=:speed)
hmξ = heatmap!(axξ, ξn, colorrange=ξ_limits, colormap=:balance)
hmθ = heatmap!(axθ, θn, colorrange=θ_limits, colormap=:thermal)
hmq = heatmap!(axq, qᵗn, colorrange=(0, qᵗ_max), colormap=Reverse(:Purples_4))
hmT = heatmap!(axT, Tn, colorrange=T_limits)
hmqˡ = heatmap!(axqˡ, qˡn, colorrange=(0, qˡ_max), colormap=Reverse(:Blues_4))Heatmap{Tuple{Vector{Float64}, Vector{Float64}, Matrix{Float32}}}Plot the surface fluxes
lines!(axτ, τˣn, color=:black, linewidth=2)
lines!(ax𝒬, 𝒬ᵀn, color=:firebrick, linewidth=2, label="sensible")
lines!(ax𝒬, 𝒬ᵛn, color=:blue, linewidth=2, label="latent")
lines!(ax𝒬, Σ𝒬n, color=:green, linewidth=4, label="total")
Legend(fig[4, 3], ax𝒬)Legend()Add zero lines, fix axis limits, and add colorbars.
for ax in (axτ, ax𝒬)
lines!(ax, [-grid.Lx/2, grid.Lx/2], [0, 0], color=:grey, linestyle=:dash)
end
for ax in (axs, axξ, axθ, axq, axT, axqˡ, axτ, ax𝒬)
xlims!(ax, -grid.Lx/2, grid.Lx/2)
end
ylims!(axτ, -τˣ_max, τˣ_max)
ylims!(ax𝒬, 𝒬_min, 𝒬_max)
Colorbar(fig[1, 0], hms, label="√(u² + w²) (m/s)", flipaxis=false)
Colorbar(fig[1, 3], hmξ, label="∂u/∂z - ∂w/∂x (s⁻¹)")
Colorbar(fig[2, 0], hmθ, label="θ (K)", flipaxis=false)
Colorbar(fig[2, 3], hmq, label="qᵗ (kg/kg)")
Colorbar(fig[3, 0], hmT, label="T (K)", flipaxis=false)
Colorbar(fig[3, 3], hmqˡ, label="qˡ (kg/kg)")Colorbar()Now we are ready to make a cool animation.
CairoMakie.record(fig, "prescribed_sea_surface_temperature.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.2
Commit ca9b6662be4 (2025-11-20 16:25 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.2
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.2.1 `.`
[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.10
[9e8cae18] Oceananigans v0.103.1
[a01a1ee8] RRTMGP v0.21.6
[b77e0a4c] InteractiveUtils v1.11.0
[44cfe95a] Pkg v1.12.0
[9a3f8284] Random v1.11.0
This page was generated using Literate.jl.