Skip to content

Commit 0e88889

Browse files
committed
Account for ERA5/ClimaAtmos topo differences when initializing pressure. Account for lapse rate in mslp sea level reduction (diagnostics).
1 parent f7b702b commit 0e88889

File tree

3 files changed

+148
-18
lines changed

3 files changed

+148
-18
lines changed

src/diagnostics/core_diagnostics.jl

Lines changed: 21 additions & 13 deletions
Original file line numberDiff line numberDiff line change
@@ -1647,24 +1647,32 @@ add_diagnostic_variable!(
16471647
function compute_mslp!(out, state, cache, time)
16481648
thermo_params = CAP.thermodynamics_params(cache.params)
16491649
g = TD.Parameters.grav(thermo_params)
1650-
R_m_surf = Fields.level(
1651-
lazy.(TD.gas_constant_air.(thermo_params, cache.precomputed.ᶜts)),
1652-
1,
1653-
)
1650+
ts_level = Fields.level(cache.precomputed.ᶜts, 1)
1651+
R_m_surf = @. lazy(TD.gas_constant_air(thermo_params, ts_level))
16541652

1655-
# get pressure, temperature, and height at the lowest atmospheric level
16561653
p_level = Fields.level(cache.precomputed.ᶜp, 1)
1657-
t_level = Fields.level(
1658-
lazy.(TD.air_temperature.(thermo_params, cache.precomputed.ᶜts)),
1659-
1,
1660-
)
1654+
t_level = @. lazy(TD.air_temperature(thermo_params, ts_level))
16611655
z_level = Fields.level(Fields.coordinate_field(state.c.ρ).z, 1)
16621656

1663-
# compute sea level pressure using the hypsometric equation
1657+
# Reduce to mean sea level using hypsometric formulation with lapse rate adjustment
1658+
# Using constant lapse rate Γ = 6.5 K/km, with virtual temperature
1659+
# represented via R_m_surf. This reduces biases over
1660+
# very cold or very warm high-topography regions.
1661+
FT = Spaces.undertype(Fields.axes(state.c.ρ))
1662+
Γ = FT(6.5e-3) # K m^-1
1663+
1664+
# p_msl = p_z0 * [1 + Γ * z / T_z0]^( g / (R_m Γ))
1665+
# where:
1666+
# - p_z0 pressure at the lowest model level
1667+
# - T_z0 air temperature at the lowest model level
1668+
# - R_m moist-air gas constant at the surface (R_m_surf), which
1669+
# accounts for virtual-temperature effects in the exponent
1670+
# - Γ constant lapse rate (6.5 K/km here)
1671+
16641672
if isnothing(out)
1665-
return @. p_level * exp(g * z_level / (R_m_surf * t_level))
1673+
return p_level .* (1 .+ Γ .* z_level ./ t_level) .^ (g / Γ ./ R_m_surf)
16661674
else
1667-
@. out = p_level * exp(g * z_level / (R_m_surf * t_level))
1675+
out .= p_level .* (1 .+ Γ .* z_level ./ t_level) .^ (g / Γ ./ R_m_surf)
16681676
end
16691677
end
16701678

@@ -1673,7 +1681,7 @@ add_diagnostic_variable!(
16731681
long_name = "Mean Sea Level Pressure",
16741682
standard_name = "mean_sea_level_pressure",
16751683
units = "Pa",
1676-
comments = "Mean sea level pressure computed from the hypsometric equation",
1684+
comments = "Mean sea level pressure computed using a lapse-rate-dependent hypsometric reduction (ERA-style; Γ=6.5 K/km with virtual temperature via moist gas constant).",
16771685
compute! = compute_mslp!,
16781686
)
16791687

src/initial_conditions/initial_conditions.jl

Lines changed: 106 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -386,6 +386,92 @@ function overwrite_initial_conditions!(
386386
)
387387
end
388388

389+
"""
390+
correct_surface_pressure_for_topography!(
391+
p_sfc,
392+
file_path,
393+
face_space,
394+
Y,
395+
ᶜT,
396+
ᶜq_tot,
397+
thermo_params,
398+
regridder_kwargs;
399+
surface_altitude_var = "z_sfc",
400+
)
401+
402+
Adjusts the surface pressure field `p_sfc` to account for mismatches between
403+
ERA5 (file) surface altitude and the model orography when specifying pressure.
404+
405+
Δz = z_model_surface - z_sfc
406+
407+
and applies a hydrostatic correction at the surface using the local moist gas
408+
constant and temperature at the surface:
409+
410+
p_sfc .= p_sfc .* exp.(-Δz * g ./ (R_m_sfc .* T_sfc))
411+
412+
where:
413+
- `g` is gravitational acceleration from `thermo_params`
414+
- `R_m_sfc` is the moist-air gas constant evaluated from `ᶜq_tot` at the surface
415+
- `T_sfc` is the air temperature from `ᶜT` at the surface
416+
417+
Returns `true` if the correction is applied; returns `false` if the surface
418+
altitude field cannot be loaded.
419+
420+
Arguments
421+
- `p_sfc`: face field of surface pressure to be corrected (modified in-place)
422+
- `file_path`: path to the ERA5-derived initialization NetCDF file
423+
- `face_space`: face space of the model grid (for reading/regridding)
424+
- `Y`: prognostic state, used to obtain model surface height
425+
- `ᶜT`: center field of temperature
426+
- `ᶜq_tot`: center field of total specific humidity
427+
- `thermo_params`: thermodynamics parameter set
428+
- `regridder_kwargs`: keyword arguments forwarded to the regridder
429+
- `surface_altitude_var`: variable name for surface altitude (default `"z_sfc"`)
430+
"""
431+
function correct_surface_pressure_for_topography!(
432+
p_sfc,
433+
file_path,
434+
face_space,
435+
Y,
436+
ᶜT,
437+
ᶜq_tot,
438+
thermo_params,
439+
regridder_kwargs;
440+
surface_altitude_var = "z_sfc",
441+
)
442+
ᶠz_surface = Fields.level(
443+
SpaceVaryingInputs.SpaceVaryingInput(
444+
file_path,
445+
surface_altitude_var,
446+
face_space,
447+
regridder_kwargs = regridder_kwargs,
448+
),
449+
Fields.half,
450+
)
451+
452+
if ᶠz_surface === nothing
453+
return false
454+
end
455+
456+
FT = eltype(thermo_params)
457+
grav = thermo_params.grav
458+
459+
ᶠz_model_surface = Fields.level(Fields.coordinate_field(Y.f).z, Fields.half)
460+
ᶠΔz = zeros(face_space)
461+
@. ᶠΔz = ᶠz_model_surface - ᶠz_surface
462+
463+
ᶠR_m = ᶠinterp.(TD.gas_constant_air.(thermo_params, TD.PhasePartition.(ᶜq_tot)))
464+
ᶠR_m_sfc = Fields.level(ᶠR_m, Fields.half)
465+
466+
ᶠT = ᶠinterp.(ᶜT)
467+
ᶠT_sfc = Fields.level(ᶠT, Fields.half)
468+
469+
@. p_sfc = p_sfc * exp(FT(-1) * ᶠΔz * grav / (ᶠR_m_sfc * ᶠT_sfc))
470+
471+
@info "Adjusted surface pressure to account for ERA5/model surface-height differences."
472+
return true
473+
end
474+
389475
# WeatherModel function using the shared implementation
390476
function overwrite_initial_conditions!(
391477
initial_condition::WeatherModel,
@@ -436,6 +522,26 @@ function overwrite_initial_conditions!(
436522
center_space,
437523
regridder_kwargs = regridder_kwargs,
438524
)
525+
# Apply hydrostatic surface-pressure correction only if surface altitude is available
526+
surface_altitude_var = "z_sfc"
527+
has_surface_altitude = NC.NCDataset(file_path) do ds
528+
haskey(ds, surface_altitude_var)
529+
end
530+
if has_surface_altitude
531+
correct_surface_pressure_for_topography!(
532+
p_sfc,
533+
file_path,
534+
face_space,
535+
Y,
536+
ᶜT,
537+
ᶜq_tot,
538+
thermo_params,
539+
regridder_kwargs;
540+
surface_altitude_var = surface_altitude_var,
541+
)
542+
else
543+
@warn "Skipping topographic correction because variable `$surface_altitude_var` is missing from $(file_path)."
544+
end
439545

440546
# With the known temperature (ᶜT) and moisture (ᶜq_tot) profile,
441547
# recompute the pressure levels assuming hydrostatic balance is maintained.

src/utils/weather_model.jl

Lines changed: 21 additions & 5 deletions
Original file line numberDiff line numberDiff line change
@@ -158,13 +158,29 @@ function to_z_levels(source_file, target_file, target_levels, FT)
158158
end
159159

160160
# Write 2D surface variables - extend to all levels (TODO: accept 2D variables in atmos)
161-
# Simply repeat the surface values for all levels
162-
surf_map = Dict("skt" => "skt", "sp" => "p")
161+
# Duplicate 2D surface field across all target vertical levels
162+
surf_map = Dict("skt" => "skt", "sp" => "p", "surface_geopotential" => "z_sfc")
163163
for (src_name, dst_name) in surf_map
164-
var_obj =
165-
defVar(ncout, dst_name, FT, ("lon", "lat", "z"), attrib = ncin[src_name].attrib)
164+
# Choose attributes; for z_sfc, set clean altitude attributes
165+
var_attrib = if dst_name == "z_sfc"
166+
Dict(
167+
"standard_name" => "surface_altitude",
168+
"long_name" => "surface altitude derived from ERA5",
169+
"units" => "m",
170+
"source_variable" => src_name,
171+
)
172+
else
173+
ncin[src_name].attrib
174+
end
175+
var_obj = defVar(ncout, dst_name, FT, ("lon", "lat", "z"), attrib = var_attrib)
176+
# Read first time slice and coalesce; follow same convention as sp (use [:, :, 1])
177+
data2d = FT.(coalesce.(ncin[src_name][:, :, 1], NaN))
178+
# Convert geopotential to meters if necessary
179+
if dst_name == "z_sfc"
180+
data2d .= data2d ./ grav
181+
end
166182
for k in 1:length(target_levels)
167-
var_obj[:, :, k] = FT.(ncin[src_name][:, :, 1])
183+
var_obj[:, :, k] = data2d
168184
end
169185
end
170186

0 commit comments

Comments
 (0)