Schär mountain wave with terrain-following coordinates
This example simulates the classic Schär mountain wave test case using the fully compressible equations on a terrain-following grid. The test case features a bell-shaped mountain with superimposed small-scale corrugations, generating both propagating and evanescent wave modes.
The Schär mountain wave (Schär et al. (2002)) is a stringent benchmark for terrain-following coordinates because the fine-scale terrain corrugations create steep coordinate-surface slopes that challenge the horizontal pressure gradient computation. We deform the computational grid itself using follow_terrain!, which applies the Gal-Chen and Somerville (1975) coordinate transformation to exactly capture the mountain profile at any resolution.
References
- Gal-Chen, T. and Somerville, R. C. (1975). On the use of a coordinate transformation for the solution of the Navier-Stokes equations. Journal of Computational Physics 17, 209–228.
- Klemp, J. B. (2011). A terrain-following coordinate with smoothed coordinate surfaces. Monthly Weather Review 139, 2163–2169.
- Klemp, J. B.; Skamarock, W. C. and Park, S.-H. (2015). Idealized global nonhydrostatic atmospheric test cases on a reduced-radius sphere. Journal of Advances in Modeling Earth Systems 7, 1155–1177.
- Schär, C.; Leuenberger, D.; Fuhrer, O.; Lüthi, D. and Girard, C. (2002). A new terrain-following vertical coordinate formulation for atmospheric prediction models. Monthly Weather Review 130, 2459–2480.
Physical setup
The simulation initializes an atmosphere with constant Brunt–Väisälä frequency $N$ and uniform background wind $U$. A bell-shaped mountain with superimposed oscillations triggers both propagating gravity waves (for wavenumbers below the critical value $k^*$) and evanescent waves (above $k^*$).
Mountain profile
The terrain follows the Schär mountain profile (Equation 46 by Schär et al. (2002)):
\[h(x) = h_0 \exp\left(-\frac{x^2}{a^2}\right) \cos^2\left(\frac{\pi x}{\lambda}\right)\]
where $h_0$ is the peak height, $a$ is the Gaussian half-width, and $\lambda$ is the wavelength of the terrain corrugations.
Constant stratification base state
For a reference temperature $T_0$, the background atmosphere corresponds to a constant Brunt–Väisälä frequency $N$:
\[N^2 = \frac{g^2}{c_p^d T_0}\]
The density-scale height parameter is:
\[\beta = \frac{g}{R^d T_0}\]
The potential temperature profile that maintains constant $N$ is:
\[\theta(z) = \theta_0 \exp\left(\frac{N^2 z}{g}\right)\]
Linear wave theory
For the linearized mountain wave problem, vertical wavenumber $m$ satisfies the dispersion relation (Appendix A of Klemp et al. (2015)):
\[m^2 = \frac{N^2}{U^2} - \frac{\beta^2}{4} - k^2\]
Waves propagate vertically when $m^2 > 0$ (i.e., $k < k^*$), and decay exponentially when $m^2 < 0$ (i.e., $k > k^*$), where the critical wavenumber is:
\[k^* = \sqrt{\frac{N^2}{U^2} - \frac{\beta^2}{4}}\]
Open lateral boundary conditions have not been implemented; periodic boundaries are used instead, which is not ideal for this test case.
using Breeze
using Oceananigans
using Breeze.TerrainFollowingDiscretization: SlopeInsideInterpolation
using Oceananigans.Grids: MutableVerticalDiscretization, xnodes, znode
using Oceananigans.Operators: Δzᶜᶜᶠ
using Oceananigans.Units
using Oceananigans: Face, Center
using Breeze.Thermodynamics: dry_air_gas_constant, hydrostatic_pressure
using Printf
using CairoMakie
using CUDAThermodynamic parameters
We define the base state with surface pressure $p_0 = 1000 \, {\rm hPa}$ and reference temperature $T_0 = 300 \, {\rm K}$. We also set the background wind at $U = 20 \, {\rm m/s}$:
constants = ThermodynamicConstants(Float64)
g = constants.gravitational_acceleration
cᵖᵈ = constants.dry_air.heat_capacity
Rᵈ = dry_air_gas_constant(constants)
p₀ = 100000 # Pa - surface pressure
T₀ = 300 # K - reference temperature
θ₀ = T₀ # K - reference potential temperature
U = 20 # m s⁻¹ - uniform background wind20Derived atmospheric parameters for an isothermal base state at $T_0$. Note: the standard Schär et al. (2002) parameters use $N = 0.01 \, {\rm s}^{-1}$ (set N² = 1e-4 to match that paper exactly).
N² = g^2 / (cᵖᵈ * T₀) # s⁻² - Brunt–Väisälä frequency squared (≈ 3.2e-4)
N = sqrt(N²) # s⁻¹ - Brunt–Väisälä frequency
β = g / (Rᵈ * T₀) # m⁻¹ - density scale parameter0.0001139362871088201Schär mountain parameters
The mountain profile with the parameters used by Schär et al. (2002) is:
h₀ = 250 # m - peak mountain height (use 25 m for strict linearity)
a = 5000 # m - Gaussian half-width parameter
λ = 4000 # m - terrain corrugation wavelength
K = 2π / λ # rad m⁻¹ - terrain wavenumber
hill(x, y) = h₀ * exp(-(x / a)^2) * cos(π * x / λ)^2hill (generic function with 1 method)Grid setup
The domain spans 200 km horizontally (±100 km) and 20 km vertically. We use a MutableVerticalDiscretization with uniform spacing in the computational coordinate; the terrain-following transformation deforms these surfaces to follow the mountain, so no vertical grid stretching near the surface is needed.
arch = CUDA.functional() ? GPU() : CPU()
Nx, Nz = 201, 50
L = 100kilometers
const H = 20kilometers
z_faces = MutableVerticalDiscretization(collect(range(0, H, length=Nz+1)))
grid = RectilinearGrid(arch, size = (Nx, Nz),
x = (-L/2, L/2), z = z_faces,
topology = (Periodic, Flat, Bounded))201×1×50 RectilinearGrid{Float64, Periodic, Flat, Bounded} on CUDAGPU with 3×0×3 halo
├── Periodic x ∈ [-50000.0, 50000.0) regularly spaced with Δx=497.512
├── Flat y
└── Bounded z ∈ [0.0, 20000.0] variably and mutably spaced with min(Δr)=400.0, max(Δr)=400.0Terrain
Apply the Schär mountain profile to the grid via follow_terrain!. This sets the column scaling factors $\sigma$ and surface displacement $\eta$ on the grid's MutableVerticalDiscretization and returns a TerrainMetrics object with precomputed terrain slopes.
metrics = follow_terrain!(grid, hill; pressure_gradient_stencil = SlopeInsideInterpolation())TerrainMetrics
├── topography: 201×1×1 Field{Oceananigans.Grids.Center, Oceananigans.Grids.Center, Oceananigans.Grids.Center} on Oceananigans.Grids.RectilinearGrid on CUDAGPU
├── ∂x_h: 201×1×1 Field{Oceananigans.Grids.Face, Oceananigans.Grids.Center, Oceananigans.Grids.Center} on Oceananigans.Grids.RectilinearGrid on CUDAGPU
├── ∂y_h: 201×1×1 Field{Oceananigans.Grids.Center, Oceananigans.Grids.Face, Oceananigans.Grids.Center} on Oceananigans.Grids.RectilinearGrid on CUDAGPU
├── z_top: 20000.0
└── pressure_gradient_stencil: Breeze.TerrainFollowingDiscretization.SlopeInsideInterpolationPlot: Terrain-following coordinate surfaces
Visualize how the computational grid deforms to follow the corrugated terrain. The terrain-following coordinate exactly captures the mountain profile regardless of vertical resolution.
x_grid = xnodes(grid, Center())
fig = Figure(size=(900, 400))
ax = Axis(fig[1, 1],
xlabel = "x (m)",
ylabel = "Height (m)",
title = "Terrain-following grid near the Schär mountain")
# Plot coordinate surfaces at selected vertical levels
CUDA.@allowscalar for k in [1, 3, 5, 8, 12, 17, 25]
z_surface = [znode(i, 1, k, grid, Center(), Center(), Face()) for i in 1:Nx]
lines!(ax, x_grid, z_surface, color = :gray, linewidth = 0.5)
end
# The bottom coordinate surface (terrain)
z_bottom = CUDA.@allowscalar [znode(i, 1, 1, grid, Center(), Center(), Face()) for i in 1:Nx]
lines!(ax, x_grid, z_bottom, color = :brown, linewidth = 2, label = "Terrain surface")
band!(ax, x_grid, zero(z_bottom), z_bottom, color = (:brown, 0.2))
# Analytical terrain for comparison
x_fine = range(-L/6, L/6, length=2000)
lines!(ax, collect(x_fine), [hill(x, 0) for x in x_fine],
color = :black, linestyle = :dash, linewidth = 1, label = "Analytical h(x)")
axislegend(ax, position = :rt)
xlims!(ax, -L/6, L/6)
ylims!(ax, -100, 3500)
figPotential temperature profile
The potential temperature profile that maintains constant Brunt–Väisälä frequency is:
\[\theta(z) = \theta_0 \exp\left(\frac{N^2 z}{g}\right)\]
θ_of_z(z) = θ₀ * exp(N² * z / g)θ_of_z (generic function with 1 method)Rayleigh damping layer
A sponge layer near the domain top absorbs upward-propagating waves and prevents spurious reflections from the rigid lid. The Gaussian mask is concentrated in the upper quarter of the domain and decays to effectively zero below mid-level.
const sponge_width = H / 4
sponge_mask(x, z) = exp(-(z - H)^2 / sponge_width^2)
ρw_sponge = Relaxation(rate=1/10, mask=sponge_mask)Relaxation{Float64, typeof(Main.var"##277".sponge_mask), typeof(Oceananigans.Forcings.zerofunction)}
├── rate: 0.1
├── mask: sponge_mask (generic function with 1 method)
└── target: 0Model construction
Build a compressible model with explicit time-stepping, 5th-order WENO advection, and terrain corrections. Passing terrain_metrics to CompressibleDynamics activates the terrain-following physics: contravariant vertical velocity, corrected pressure gradient, and terrain-aware divergence. The reference_potential_temperature enables a perturbation pressure approach for the horizontal pressure gradient that reduces the truncation error inherent in terrain-following coordinates (Klemp (2011)).
pˢᵗ = 1e5 # Pa - standard pressure
κ = Rᵈ / cᵖᵈ
dynamics = CompressibleDynamics(ExplicitTimeStepping();
terrain_metrics = metrics,
surface_pressure = p₀,
reference_potential_temperature = θ_of_z)
model = AtmosphereModel(grid; dynamics, advection=WENO(order=5), forcing=(; ρw=ρw_sponge))AtmosphereModel{GPU, RectilinearGrid}(time = 0 seconds, iteration = 0)
├── grid: 201×1×50 RectilinearGrid{Float64, Periodic, Flat, Bounded} on CUDAGPU with 3×0×3 halo
├── dynamics: CompressibleDynamics{ExplicitTimeStepping}
├── formulation: LiquidIcePotentialTemperatureFormulation
├── thermodynamic_constants: ThermodynamicConstants{Float64}
├── timestepper: SSPRungeKutta3
├── advection scheme:
│ ├── momentum: WENO{3, Float64, Oceananigans.Utils.BackendOptimizedDivision}(order=5)
│ ├── ρθ: WENO{3, Float64, Oceananigans.Utils.BackendOptimizedDivision}(order=5)
│ └── ρqᵛ: WENO{3, Float64, Oceananigans.Utils.BackendOptimizedDivision}(order=5)
├── forcing: ρw=>ContinuousForcing
├── tracers: ()
├── coriolis: Nothing
└── microphysics: NothingInitial conditions
We initialize the atmosphere in discrete hydrostatic balance using Exner function integration. This is essential for compressible models on terrain-following grids: the equation of state alone does not produce a pressure field in discrete hydrostatic balance. Because we passed reference_potential_temperature to CompressibleDynamics, the model has already computed this discrete reference state for us in terrain_reference_density!
ρ_ref = model.dynamics.terrain_reference_density
ρ = model.dynamics.density
ρθ = model.formulation.potential_temperature_density
set!(ρ, ρ_ref)
set!(ρθ, (x, z) -> θ_of_z(z))
parent(ρθ) .*= parent(ρ)
set!(model, u=U)Simulation
Acoustic waves require a CFL condition based on the sound speed, which gives a much smaller time step than the anelastic equations. We use the potential temperature at the domain top for a conservative estimate, since θ (and thus the sound speed) increases with height.
γ = cᵖᵈ / (cᵖᵈ - Rᵈ)
θ_top = θ_of_z(H)
ℂᵃᶜ = sqrt(γ * Rᵈ * θ_top)
Δx = L / Nx
Δz = H / Nz
Δt = 0.4 * min(Δx, Δz) / (ℂᵃᶜ + U)
stop_time = 2hours
simulation = Simulation(model; Δt, stop_time)
Oceananigans.Diagnostics.erroring_NaNChecker!(simulation)Progress callback to monitor simulation health:
wall_clock = Ref(time_ns())
function progress(sim)
elapsed = 1e-9 * (time_ns() - wall_clock[])
u, v, w = sim.model.velocities
msg = @sprintf("Iter: %d, time: %s, wall time: %s, max|u|: %.2f, max|w|: %.4f m/s",
iteration(sim), prettytime(sim), prettytime(elapsed),
maximum(abs, u), maximum(abs, w))
wall_clock[] = time_ns()
@info msg
return nothing
end
add_callback!(simulation, progress, name=:progress, IterationInterval(500))Output
Save vertical velocity and contravariant vertical velocity for post-processing:
w = model.velocities.w
Ω̃ = model.dynamics.contravariant_vertical_velocity
filename = "mountain_waves"
outputs = (; w, Ω̃)
simulation.output_writers[:jld2] = JLD2Writer(model, outputs;
filename,
schedule = TimeInterval(2minutes),
overwrite_existing = true)
run!(simulation)[ Info: Initializing simulation...
[ Info: Iter: 0, time: 0 seconds, wall time: 25.183 seconds, max|u|: 20.00, max|w|: 0.0000 m/s
[ Info: ... simulation initialization complete (32.468 seconds)
[ Info: Executing initial time step...
[ Info: ... initial time step complete (7.554 seconds).
[ Info: Iter: 500, time: 2.660 minutes, wall time: 16.813 seconds, max|u|: 22.69, max|w|: 2.3591 m/s
[ Info: Iter: 1000, time: 5.321 minutes, wall time: 2.820 seconds, max|u|: 22.61, max|w|: 2.6626 m/s
[ Info: Iter: 1500, time: 7.981 minutes, wall time: 2.774 seconds, max|u|: 22.87, max|w|: 2.5738 m/s
[ Info: Iter: 2000, time: 10.639 minutes, wall time: 2.826 seconds, max|u|: 22.88, max|w|: 2.5855 m/s
[ Info: Iter: 2500, time: 13.300 minutes, wall time: 2.826 seconds, max|u|: 22.89, max|w|: 2.5714 m/s
[ Info: Iter: 3000, time: 15.960 minutes, wall time: 2.840 seconds, max|u|: 23.07, max|w|: 2.5887 m/s
[ Info: Iter: 3500, time: 18.618 minutes, wall time: 2.834 seconds, max|u|: 22.99, max|w|: 2.6091 m/s
[ Info: Iter: 4000, time: 21.278 minutes, wall time: 2.817 seconds, max|u|: 22.96, max|w|: 2.5879 m/s
[ Info: Iter: 4500, time: 23.939 minutes, wall time: 2.840 seconds, max|u|: 22.99, max|w|: 2.5981 m/s
[ Info: Iter: 5000, time: 26.597 minutes, wall time: 3.514 seconds, max|u|: 23.04, max|w|: 2.5810 m/s
[ Info: Iter: 5500, time: 29.257 minutes, wall time: 2.830 seconds, max|u|: 23.07, max|w|: 2.5896 m/s
[ Info: Iter: 6000, time: 31.917 minutes, wall time: 2.804 seconds, max|u|: 22.76, max|w|: 2.5575 m/s
[ Info: Iter: 6500, time: 34.575 minutes, wall time: 2.836 seconds, max|u|: 23.16, max|w|: 2.6089 m/s
[ Info: Iter: 7000, time: 37.236 minutes, wall time: 2.831 seconds, max|u|: 22.96, max|w|: 2.6003 m/s
[ Info: Iter: 7500, time: 39.896 minutes, wall time: 2.809 seconds, max|u|: 23.17, max|w|: 2.5767 m/s
[ Info: Iter: 8000, time: 42.554 minutes, wall time: 2.829 seconds, max|u|: 23.02, max|w|: 2.5800 m/s
[ Info: Iter: 8500, time: 45.214 minutes, wall time: 2.812 seconds, max|u|: 22.99, max|w|: 2.5589 m/s
[ Info: Iter: 9000, time: 47.875 minutes, wall time: 2.829 seconds, max|u|: 22.99, max|w|: 2.5768 m/s
[ Info: Iter: 9500, time: 50.533 minutes, wall time: 2.820 seconds, max|u|: 23.02, max|w|: 2.6039 m/s
[ Info: Iter: 10000, time: 53.193 minutes, wall time: 2.824 seconds, max|u|: 22.97, max|w|: 2.5815 m/s
[ Info: Iter: 10500, time: 55.854 minutes, wall time: 2.832 seconds, max|u|: 23.19, max|w|: 2.5997 m/s
[ Info: Iter: 11000, time: 58.511 minutes, wall time: 2.844 seconds, max|u|: 23.03, max|w|: 2.5694 m/s
[ Info: Iter: 11500, time: 1.020 hours, wall time: 2.831 seconds, max|u|: 23.22, max|w|: 2.5939 m/s
[ Info: Iter: 12000, time: 1.064 hours, wall time: 2.825 seconds, max|u|: 22.56, max|w|: 2.5525 m/s
[ Info: Iter: 12500, time: 1.108 hours, wall time: 2.842 seconds, max|u|: 23.10, max|w|: 2.5923 m/s
[ Info: Iter: 13000, time: 1.153 hours, wall time: 3.527 seconds, max|u|: 23.06, max|w|: 2.5958 m/s
[ Info: Iter: 13500, time: 1.197 hours, wall time: 2.856 seconds, max|u|: 23.04, max|w|: 2.5769 m/s
[ Info: Iter: 14000, time: 1.241 hours, wall time: 2.853 seconds, max|u|: 23.01, max|w|: 2.5713 m/s
[ Info: Iter: 14500, time: 1.285 hours, wall time: 2.857 seconds, max|u|: 22.79, max|w|: 2.5577 m/s
[ Info: Iter: 15000, time: 1.330 hours, wall time: 2.828 seconds, max|u|: 22.90, max|w|: 2.5749 m/s
[ Info: Iter: 15500, time: 1.374 hours, wall time: 2.846 seconds, max|u|: 23.06, max|w|: 2.6160 m/s
[ Info: Iter: 16000, time: 1.418 hours, wall time: 2.849 seconds, max|u|: 22.98, max|w|: 2.5343 m/s
[ Info: Iter: 16500, time: 1.463 hours, wall time: 2.840 seconds, max|u|: 22.90, max|w|: 2.5923 m/s
[ Info: Iter: 17000, time: 1.507 hours, wall time: 2.832 seconds, max|u|: 22.95, max|w|: 2.5629 m/s
[ Info: Iter: 17500, time: 1.551 hours, wall time: 2.860 seconds, max|u|: 22.89, max|w|: 2.5736 m/s
[ Info: Iter: 18000, time: 1.596 hours, wall time: 2.846 seconds, max|u|: 22.65, max|w|: 2.5492 m/s
[ Info: Iter: 18500, time: 1.640 hours, wall time: 2.852 seconds, max|u|: 22.99, max|w|: 2.5671 m/s
[ Info: Iter: 19000, time: 1.684 hours, wall time: 2.837 seconds, max|u|: 23.04, max|w|: 2.5892 m/s
[ Info: Iter: 19500, time: 1.729 hours, wall time: 2.847 seconds, max|u|: 22.98, max|w|: 2.5729 m/s
[ Info: Iter: 20000, time: 1.773 hours, wall time: 2.853 seconds, max|u|: 22.99, max|w|: 2.5789 m/s
[ Info: Iter: 20500, time: 1.817 hours, wall time: 2.835 seconds, max|u|: 22.79, max|w|: 2.5503 m/s
[ Info: Iter: 21000, time: 1.862 hours, wall time: 2.824 seconds, max|u|: 23.14, max|w|: 2.5833 m/s
[ Info: Iter: 21500, time: 1.906 hours, wall time: 3.373 seconds, max|u|: 23.13, max|w|: 2.6020 m/s
[ Info: Iter: 22000, time: 1.950 hours, wall time: 2.708 seconds, max|u|: 23.17, max|w|: 2.5915 m/s
[ Info: Iter: 22500, time: 1.995 hours, wall time: 2.808 seconds, max|u|: 23.02, max|w|: 2.5792 m/s
[ Info: Simulation is stopping after running for 2.835 minutes.
[ Info: Simulation time 2 hours equals or exceeds stop time 2 hours.
Analytical solution
The linear analytical solution for mountain waves provides a validation benchmark. Following Appendix A of Klemp et al. (2015), the vertical velocity field is computed via Fourier integration over wavenumber space.
Fourier transform of terrain
The Fourier transform of the Schär mountain profile (Equation A8):
\[\hat{h}(k) = \frac{\sqrt{\pi} h_0 a}{4} \left[ e^{-a^2(K+k)^2/4} + e^{-a^2(K-k)^2/4} + 2e^{-a^2 k^2/4} \right]\]
ĥ(k) = sqrt(π) * h₀ * a / 4 * (exp(-a^2 * (K + k)^2 / 4) +
exp(-a^2 * (K - k)^2 / 4) +
2exp(-a^2 * k^2 / 4))ĥ (generic function with 1 method)Dispersion relation
Vertical wavenumber squared (Equation A5) and critical wavenumber (Equation A11) by Klemp et al. (2015):
m²(k) = (N² / U^2 - β^2 / 4) - k^2
k★ = sqrt(N² / U^2 - β^2 / 4)0.00089147756144519Linear vertical velocity
Compute the analytical vertical velocity $w(x, z)$ via Equation A10 by Klemp et al. (2015):
\[w(x, z) = -\frac{U}{\pi} e^{\beta z/2} \left[ \int_0^{k^*} k \hat{h}(k) \sin(m z + k x) \, \mathrm{d}k + \int_{k^*}^{\infty} k \hat{h}(k) e^{-|m| z} \sin(k x) \, \mathrm{d}k \right]\]
Above, the first integral represents propagating waves and the second represents evanescent waves.
"""
w_linear(x, z; nk=100)
Compute the 2-D linear vertical velocity `w(x,z)` from the analytical solution
(Appendix A, Equation A10 of Klemp et al., 2015).
"""
function w_linear(x, z; nk=100)
k = range(0, 10k★; length=nk)
m_abs = @. sqrt(abs(m²(k)))
integrand = @. k * ĥ(k) * ifelse(m²(k) ≥ 0,
sin(m_abs * z + k * x),
exp(-m_abs * z) * sin(k * x))
# Numerical integration using trapezoidal rule:
Δk = step(k)
integral = Δk * (sum(integrand) - (first(integrand) + last(integrand)) / 2)
return -(U / π) * exp(β * z / 2) * integral
endResults: Comparison with analytical solution
We compare the simulated vertical velocity field at 10 minutes with the linear analytical solution. The terrain-following grid exactly represents the corrugated mountain profile, avoiding the staircase artifacts of the immersed boundary method. The heatmaps display fields in physical $(x, z)$ coordinates, with the deformed grid automatically mapped to the true terrain geometry.
First, we compute analytical vertical velocity $w$ on the same grid as the simulation.
w_analytical = Field{Center, Nothing, Face}(grid)
set!(w_analytical, w_linear)
fig = Figure(size=(900, 800), fontsize=14)
ax1 = Axis(fig[1, 1], ylabel = "z (m)", title = "Simulated w at t = 2 hours")
ax2 = Axis(fig[2, 1], xlabel = "x (m)", ylabel = "z (m)", title = "Linear Analytical w")
hidexdecorations!(ax1, grid = false)
hm1 = heatmap!(ax1, w, colormap = :balance, colorrange = (-1, 1))
hm2 = heatmap!(ax2, w_analytical, colormap = :balance, colorrange = (-1, 1))
Colorbar(fig[1:2, 2], hm1, label = "w (m s⁻¹)")
for ax in (ax1, ax2)
ax.limits = ((-30e3, 30e3), (0, 10e3))
end
figJulia version and environment information
This example was executed with the following version of Julia:
using InteractiveUtils: versioninfo
versioninfo()Julia Version 1.12.6
Commit 15346901f00 (2026-04-09 19:20 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.6
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.7 `.`
⌃ [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.107.3
[a01a1ee8] RRTMGP v0.21.7
[b77e0a4c] InteractiveUtils v1.11.0
[44cfe95a] Pkg v1.12.1
[9a3f8284] Random v1.11.0
Info Packages marked with ⌃ and ⌅ have new versions available. Those with ⌃ may be upgradable, but those with ⌅ are restricted by compatibility constraints from upgrading. To see why use `status --outdated`
This page was generated using Literate.jl.