Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Don't output diagnostics when simulation is less than 1h #3337

Merged
merged 1 commit into from
Sep 27, 2024
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
10 changes: 5 additions & 5 deletions perf/flame.jl
Original file line number Diff line number Diff line change
Expand Up @@ -38,13 +38,13 @@ ProfileCanvas.html_file(joinpath(output_dir, "flame.html"), results)
#####

allocs_limit = Dict()
allocs_limit["flame_perf_target"] = 1_674_376
allocs_limit["flame_perf_target_tracers"] = 1_820_752
allocs_limit["flame_perf_target"] = 33_096
allocs_limit["flame_perf_target_tracers"] = 49_264
allocs_limit["flame_perf_diagnostics"] = 12_301_560
allocs_limit["flame_perf_target_diagnostic_edmfx"] = 2_555_280
allocs_limit["flame_perf_target_diagnostic_edmfx"] = 47_744
allocs_limit["flame_perf_target_frierson"] = 1_849_976
allocs_limit["flame_perf_target_threaded"] = 2_306_856
allocs_limit["flame_perf_target_callbacks"] = 1_879_184
allocs_limit["flame_perf_target_threaded"] = 883_784
allocs_limit["flame_perf_target_callbacks"] = 230_192
allocs_limit["flame_perf_gw"] = 882_938_744
allocs_limit["flame_perf_target_prognostic_edmfx_aquaplanet"] = 2_490_248
allocs_limit["flame_gpu_implicit_barowave_moist"] = 2_099_536
Expand Down
4 changes: 2 additions & 2 deletions src/callbacks/get_callbacks.jl
Original file line number Diff line number Diff line change
@@ -1,4 +1,4 @@
function get_diagnostics(parsed_args, atmos_model, Y, p, dt)
function get_diagnostics(parsed_args, atmos_model, Y, p, dt, t_start)

FT = Spaces.undertype(axes(Y.c))

Expand Down Expand Up @@ -150,7 +150,7 @@ function get_diagnostics(parsed_args, atmos_model, Y, p, dt)
diagnostics = [
CAD.default_diagnostics(
atmos_model,
time_to_seconds(parsed_args["t_end"]),
time_to_seconds(parsed_args["t_end"]) - t_start,
p.start_date;
output_writer = netcdf_writer,
)...,
Expand Down
85 changes: 44 additions & 41 deletions src/diagnostics/default_diagnostics.jl
Original file line number Diff line number Diff line change
Expand Up @@ -5,11 +5,11 @@
# level interfaces, add them here. Feel free to include extra files.

"""
default_diagnostics(model, t_end, start_date; output_writer)
default_diagnostics(model, duration, start_date; output_writer)

Return a list of `ScheduledDiagnostic`s associated with the given `model` that use
`output_write` to write to disk. `t_end` is the expected simulation end time and it is used
to choose the most reasonable output frequency.
`output_write` to write to disk. `duration` is the expected duration of the simulation and
it is used to choose the most reasonable output frequency.

`start_date` is the date that we assign at the start of the simulation.
We convert time to date as
Expand All @@ -20,14 +20,14 @@ current_date = start_date + integrator.t

The logic is as follows:

If `t_end < 1 day` take hourly means,
if `t_end < 30 days` take daily means,
if `t_end < 90 days` take means over ten days,
If `t_end >= 90 year` take monthly means.
If `duration < 1 day` take hourly means,
if `duration < 30 days` take daily means,
if `duration < 90 days` take means over ten days,
If `duration >= 90 year` take monthly means.
"""
function default_diagnostics(
model::AtmosModel,
t_end,
duration,
start_date::DateTime;
output_writer,
)
Expand All @@ -38,7 +38,7 @@ function default_diagnostics(
x ->
default_diagnostics(
getfield(model, x),
t_end,
duration,
start_date;
output_writer,
) != [],
Expand All @@ -48,11 +48,11 @@ function default_diagnostics(
# We use a map because we want to ensure that diagnostics is a well defined type, not
# Any. This reduces latency.
return vcat(
core_default_diagnostics(output_writer, t_end, start_date),
core_default_diagnostics(output_writer, duration, start_date),
map(non_empty_fields) do field
default_diagnostics(
getfield(model, field),
t_end,
duration,
start_date;
output_writer,
)
Expand All @@ -65,7 +65,7 @@ end
# that all the default_diagnostics return the same type). This is used by
# default_diagnostics(model::AtmosModel; output_writer), so that we can ignore defaults for
# submodels that have no given defaults.
default_diagnostics(submodel, t_end, start_date; output_writer) = []
default_diagnostics(submodel, duration, start_date; output_writer) = []

"""
produce_common_diagnostic_function(period, reduction)
Expand Down Expand Up @@ -101,25 +101,28 @@ end
include("standard_diagnostic_frequencies.jl")

"""
frequency_averages(t_end::Real)
frequency_averages(duration::Real)

Return the correct averaging function depending on the total simulation time.

If `t_end < 1 day` take hourly means,
if `t_end < 30 days` take daily means,
if `t_end < 90 days` take means over ten days,
If `t_end >= 90 year` take monthly means.
If `duration < 1 hour` do nothing,
If `duration < 1 day` take hourly means,
if `duration < 30 days` take daily means,
if `duration < 90 days` take means over ten days,
If `duration >= 90 year` take monthly means.
"""
function frequency_averages(t_end)
FT = eltype(t_end)
if t_end >= 90 * 86400
function frequency_averages(duration)
FT = eltype(duration)
if duration >= 90 * 86400
return (args...; kwargs...) -> monthly_averages(FT, args...; kwargs...)
elseif t_end >= 30 * 86400
elseif duration >= 30 * 86400
return (args...; kwargs...) -> tendaily_averages(FT, args...; kwargs...)
elseif t_end >= 86400
elseif duration >= 86400
return (args...; kwargs...) -> daily_averages(FT, args...; kwargs...)
else
elseif duration >= 3600
return (args...; kwargs...) -> hourly_averages(FT, args...; kwargs...)
else
return (args...; kwargs...) -> ()
end
end

Expand All @@ -128,7 +131,7 @@ end
########
# Core #
########
function core_default_diagnostics(output_writer, t_end, start_date)
function core_default_diagnostics(output_writer, duration, start_date)
core_diagnostics = [
"ts",
"ta",
Expand All @@ -145,16 +148,16 @@ function core_default_diagnostics(output_writer, t_end, start_date)
"hfes",
]

average_func = frequency_averages(t_end)
FT = eltype(t_end)
average_func = frequency_averages(duration)
FT = eltype(duration)

if t_end >= 90 * 86400
if duration >= 90 * 86400
min_func = (args...; kwargs...) -> monthly_min(FT, args...; kwargs...)
max_func = (args...; kwargs...) -> monthly_max(FT, args...; kwargs...)
elseif t_end >= 30 * 86400
elseif duration >= 30 * 86400
min_func = (args...; kwargs...) -> tendaily_min(FT, args...; kwargs...)
max_func = (args...; kwargs...) -> tendaily_max(FT, args...; kwargs...)
elseif t_end >= 86400
elseif duration >= 86400
min_func = (args...; kwargs...) -> daily_min(FT, args...; kwargs...)
max_func = (args...; kwargs...) -> daily_max(FT, args...; kwargs...)
else
Expand Down Expand Up @@ -184,7 +187,7 @@ end
##################
function default_diagnostics(
::T,
t_end,
duration,
start_date;
output_writer,
) where {T <: Union{EquilMoistModel, NonEquilMoistModel}}
Expand All @@ -202,7 +205,7 @@ function default_diagnostics(
"clwvi",
"clivi",
]
average_func = frequency_averages(t_end)
average_func = frequency_averages(duration)
return [average_func(moist_diagnostics...; output_writer, start_date)...]
end

Expand All @@ -211,13 +214,13 @@ end
#######################
function default_diagnostics(
::Microphysics1Moment,
t_end,
duration,
start_date;
output_writer,
)
precip_diagnostics = ["husra", "hussn"]

average_func = frequency_averages(t_end)
average_func = frequency_averages(duration)

return [average_func(precip_diagnostics...; output_writer, start_date)...]
end
Expand All @@ -227,7 +230,7 @@ end
##################
function default_diagnostics(
::RRTMGPI.AbstractRRTMGPMode,
t_end,
duration,
start_date;
output_writer,
)
Expand All @@ -245,15 +248,15 @@ function default_diagnostics(
"rlus",
]

average_func = frequency_averages(t_end)
average_func = frequency_averages(duration)

return [average_func(rad_diagnostics...; output_writer, start_date)...]
end


function default_diagnostics(
::RRTMGPI.AllSkyRadiationWithClearSkyDiagnostics,
t_end,
duration,
start_date;
output_writer,
)
Expand Down Expand Up @@ -282,7 +285,7 @@ function default_diagnostics(
"rlutcs",
]

average_func = frequency_averages(t_end)
average_func = frequency_averages(duration)

return [
average_func(rad_diagnostics...; output_writer, start_date)...,
Expand All @@ -295,7 +298,7 @@ end
##################
function default_diagnostics(
::PrognosticEDMFX,
t_end,
duration,
start_date;
output_writer,
)
Expand Down Expand Up @@ -326,7 +329,7 @@ function default_diagnostics(
"lmix",
]

average_func = frequency_averages(t_end)
average_func = frequency_averages(duration)

return [
average_func(edmfx_draft_diagnostics...; output_writer, start_date)...,
Expand All @@ -337,7 +340,7 @@ end

function default_diagnostics(
::DiagnosticEDMFX,
t_end,
duration,
start_date;
output_writer,
)
Expand All @@ -355,7 +358,7 @@ function default_diagnostics(
]
diagnostic_edmfx_env_diagnostics = ["waen", "tke", "lmix"]

average_func = frequency_averages(t_end)
average_func = frequency_averages(duration)

return [
average_func(
Expand Down
1 change: 1 addition & 0 deletions src/solver/type_getters.jl
Original file line number Diff line number Diff line change
Expand Up @@ -754,6 +754,7 @@ function get_simulation(config::AtmosConfig)
Y,
p,
sim_info.dt,
t_start,
)
end
@info "initializing diagnostics: $s"
Expand Down
Loading