Global bucket run

The code sets up and runs the bucket model on a spherical domain, using ERA5 data.

First we import a lot of packages:

import ClimaComms
using ClimaCore
using ClimaUtilities
import Interpolations
import ClimaUtilities.TimeVaryingInputs:
    TimeVaryingInput, LinearInterpolation, PeriodicCalendar
ClimaComms.@import_required_backends
import ClimaTimeSteppers as CTS
import ClimaParams as CP
using ClimaLand.Bucket:
    BucketModel, BucketModelParameters, PrescribedBaregroundAlbedo
import ClimaLand
import ClimaLand.Parameters as LP
import ClimaLand.Simulations: LandSimulation, solve!
using Dates
using CairoMakie, ClimaAnalysis, GeoMakie, Poppler_jll, Printf, StatsBase
import ClimaLand.LandSimVis as LandSimVis;

Set the simulation float type, determine the context (MPI or on a single node), and device type. Create a default output directory for diagnostics.

const FT = Float64;
context = ClimaComms.context()
ClimaComms.init(context)
device = ClimaComms.device()
device_suffix = device isa ClimaComms.CPUSingleThreaded ? "cpu" : "gpu"
root_path = "bucket_longrun_$(device_suffix)"
diagnostics_outdir = joinpath(root_path, "global_diagnostics")
outdir =
    ClimaUtilities.OutputPathGenerator.generate_output_path(diagnostics_outdir);

Set timestep, startdate, stopdate:

Δt = 900.0
start_date = DateTime(2008)
stop_date = DateTime(2009);

Create the domain - this is intentionally low resolution, about 4.5 degrees x 4.5 degrees, to run quickly when making the documentation on CPU.

nelements = (20, 7)
depth = FT(3.5)
dz_tuple = FT.((1.0, 0.05))
domain =
    ClimaLand.Domains.global_domain(FT; context, nelements, depth, dz_tuple);

Parameters:

earth_param_set = LP.LandParameters(FT)
α_snow = FT(0.8)
albedo = PrescribedBaregroundAlbedo{FT}(α_snow, domain.space.surface)
bucket_parameters = BucketModelParameters(
    FT;
    albedo,
    σS_c = FT(0.2),
    W_f = FT(0.2),
    z_0m = FT(1e-3),
    z_0b = FT(1e-3),
    κ_soil = FT(1.5),
    ρc_soil = FT(2e6),
    τc = FT(float(Δt)),
);

Low-resolution forcing data from ERA5 is used here, but high-resolution should be used for production runs.

era5_ncdata_path = ClimaLand.Artifacts.era5_land_forcing_data2008_path(;
    context,
    lowres = true,
)
atmos, radiation = ClimaLand.prescribed_forcing_era5(
    era5_ncdata_path,
    domain.space.surface,
    start_date,
    earth_param_set,
    FT;
    max_wind_speed = 25.0,
    time_interpolation_method = LinearInterpolation(PeriodicCalendar()),
    regridder_type = :InterpolationsRegridder,
);

Make the model:

bucket = BucketModel(
    parameters = bucket_parameters,
    domain = domain,
    atmosphere = atmos,
    radiation = radiation,
);

Create a function which sets the initial conditions. This should have the argument structure (Y,p,t, model) in order to be used by the LandSimulation struct, below:

function set_ic!(Y, p, t, bucket)
    coords = ClimaCore.Fields.coordinate_field(Y.bucket.T)
    T_sfc_0 = 271.0
    @. Y.bucket.T = T_sfc_0 + 40 * cosd(coords.lat)^4
    Y.bucket.W .= 0.15
    Y.bucket.Ws .= 0.0
    Y.bucket.σS .= 0.0
end
set_ic! (generic function with 1 method)

Define timestepper and ODE algorithm

timestepper = CTS.RK4()
timestepper = CTS.ExplicitAlgorithm(timestepper);

Create the simulation and solve it:

simulation = LandSimulation(
    start_date,
    stop_date,
    Δt,
    bucket;
    set_ic!,
    timestepper,
    outdir,
);

solve!(simulation);
┌ Warning: overwriting diagnostic `epa` entry containing fields
│ ("nothing", "epa", "Energy per unit ground area", "energy_per_area", "J m^-2", "Vertically integrated volumetric energy per area")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `wvpa` entry containing fields
│ ("nothing", "wvpa", "Water volume per unit ground area", "water_volume_per_area", "m^3 m^-2", "Vertically integrated volumetric water per area")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `epac` entry containing fields
│ ("nothing", "epac", "Energy per unit ground area change", "energy_per_area_change", "J m^-2", "Expected change in vertically integrated volumetric energy per area")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `wvpac` entry containing fields
│ ("nothing", "wvpac", "Water volume per unit ground area change", "water_volume_per_area_change", "m^3 m^-2", "Expected change in vertically integrated volumetric water per area")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `swa` entry containing fields
│ ("nothing", "swa", "Shortwave Albedo", "sw_albedo", "", "The fraction of downwelling shortwave radiation reflected by the land surface.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `rn` entry containing fields
│ ("nothing", "rn", "Net Radiation", "net_radiation", "W m^-2", "Difference between downwelling and upwelling shortwave and longwave radiation at the land surface.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `tsfc` entry containing fields
│ ("nothing", "tsfc", "Bucket Surface Temperature", "surface_temperature", "K", "Temperature of the bucket-land surface.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `lhf` entry containing fields
│ ("nothing", "lhf", "Latent Heat Flux", "latent_heat_flux", "W m^-2", "Exchange of energy at the land-atmosphere interface due to water evaporation or sublimation.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `rae` entry containing fields
│ ("nothing", "rae", "Aerodynamic Resistance", "aerodynamic_resistance", "m s^-1", "Effiency of turbulent transport controlling the land-atmosphere exchange of sensible and latent heat.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `shf` entry containing fields
│ ("nothing", "shf", "Sensible Heat Flux", "sensible_heat_flux", "W m^-2", "Exchange of energy at the land-atmosphere interface due to temperature difference.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `vflux` entry containing fields
│ ("nothing", "vflux", "Liquid water evaporation", "vapor_flux", "m s^-1", "Flux of water from the land surface to the atmosphere. E.g., evaporation or sublimation.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `rhosfc` entry containing fields
│ ("nothing", "rhosfc", "Surface Air Density", "surface_air_density", "kg m^−3", "Density of air at the land-atmosphere interface.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `tsoil` entry containing fields
│ ("nothing", "tsoil", "Soil temperature", "soil_temperature", "K", "Soil temperature at multiple soil depth. (depth resolved)")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `wsoil` entry containing fields
│ ("nothing", "wsoil", "subsurface Water Storage", "subsurface_water_storage", "m", "Soil water content.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `wsfc` entry containing fields
│ ("nothing", "wsfc", "Surface Water Content", "surface_water_content", "m", "Water at the soil surface.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `ssfc` entry containing fields
│ ("nothing", "ssfc", "Snow Water Equivalent", "snow_water_equivalent", "m", "Snow at the soil surface, expressed in water equivalent.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `sif` entry containing fields
│ ("nothing", "sif", "Solar Induced Fluorescence", "solar_induced_fluorescence", "W m^-2", "The fluorescence of leaves induced by solar radiation at 755nm. This quantity is correlated with photosynthesis activity.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `ra` entry containing fields
│ ("nothing", "ra", "Autotrophic Respiration", "autotrophic_respiration", "mol CO2 m^-2 s^-1", "The canopy autotrophic respiration, the sum of leaves, stems and roots respiration.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `gs` entry containing fields
│ ("nothing", "gs", "Stomatal Conductance", "stomatal_conductance", "mol H2O m^-2 s^-1", "The conductance of leaves. This depends on stomatal opening. It varies with factors such as soil moisture or atmospheric water demand.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `trans` entry containing fields
│ ("nothing", "trans", "Canopy Transpiration", "canopy_transpiration", "m s^-1", "The water evaporated from the canopy due to leaf transpiration (flux of water volume, m^3 of water per m^2 of ground).")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `clhf` entry containing fields
│ ("nothing", "clhf", "Canopy Latent Heat Flux", "canopy_latent_heat_flux", "W m^-2", "The energy used for canopy transpiration.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `cshf` entry containing fields
│ ("nothing", "cshf", "Canopy Sensible Heat Flux", "canopy_sensible_heat_flux", "W m^-2", "The energy used for canopy temperature change.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `lwp` entry containing fields
│ ("nothing", "lwp", "Leaf Water Potential", "leaf_water_potential", "m", "The water potential of a leaf.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `far` entry containing fields
│ ("nothing", "far", "Root flux per ground area", "root_flux_per_ground_area", "m s^-1", "Flux of water volume per m^2 of root per second, multiplied by the area index (root area/ground area).")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `lai` entry containing fields
│ ("nothing", "lai", "Leaf area Index", "leaf_area_index", "m^2 m^-2", "The area index of leaves, expressed in surface area of leaves per surface area of ground.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `msf` entry containing fields
│ ("nothing", "msf", "Moisture Stress Factor", "moisture_stress_factor", "", "Sensitivity of plants conductance to soil water content. Unitless")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `rai` entry containing fields
│ ("nothing", "rai", "Root area Index", "root_area_index", "m^2 m^-2", "The area index of roots, expressed in surface area of roots per surface area of ground.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `sai` entry containing fields
│ ("nothing", "sai", "Stem area Index", "stem_area_index", "m^2 m^-2", "The area index of stems, expressed in surface area of stems per surface area of ground.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `gpp` entry containing fields
│ ("nothing", "gpp", "Gross Primary Productivity", "gross_primary_productivity", "mol CO2 m^-2 s^-1", "Net photosynthesis (carbon assimilation) of the canopy. This is equivalent to leaf net assimilation scaled to the canopy level.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `an` entry containing fields
│ ("nothing", "an", "Leaf Net Photosynthesis", "leaf_net_photosynthesis", "mol CO2 m^-2 s^-1", "Net photosynthesis (carbon assimilation) of a leaf, computed for example by the Farquhar model.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `rd` entry containing fields
│ ("nothing", "rd", "Leaf Respiration", "leaf_dark_respiration", "mol CO2 m^-2 s^-1", "Leaf respiration, called dark respiration because usually measured in the abscence of radiation.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `vcmax25` entry containing fields
│ ("nothing", "vcmax25", "Vcmax25", "vcmax25", "mol CO2 m^-2 s^-1", "The parameter vcmax at 25 degree celsius. Important for the Farquhar model of leaf photosynthesis.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `nir` entry containing fields
│ ("nothing", "nir", "Near Infrared Radiation", "near_infrared_radiation", "mol photons m^-2 s^-1", "The amount of near infrared radiation reaching the canopy.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `anir` entry containing fields
│ ("nothing", "anir", "Absorbed Near Infrared Radiation", "absorbed_near_infrared_radiation", "mol photons m^-2 s^-1", "The amount of near infrared radiation absorbed by the canopy.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `rnir` entry containing fields
│ ("nothing", "rnir", "Reflected Near Infrared Radiation", "reflected_near_infrared_radiation", "mol photons m^-2 s^-1", "The amount of near infrared radiation reflected by the canopy.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `tnir` entry containing fields
│ ("nothing", "tnir", "Transmitted Near Infrared Radiation", "transmitted_near_infrared_radiation", "mol photons m^-2 s^-1", "The amount of near infrared radiation transmitted by the canopy.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `par` entry containing fields
│ ("nothing", "par", "Photosynthetically Active Radiation", "photosynthetically_active_radiation", "mol photons m^-2 s^-1", "The subset of total radiation that activates photosynthesis reaching the canopy.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `apar` entry containing fields
│ ("nothing", "apar", "Absorbed Photosynthetically Active Radiation", "absorbed_photosynthetically_active_radiation", "mol photons m^-2 s^-1", "The amount of photosynthetically active radiation absorbed by the leaf. The rest if reflected or transmitted.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `rpar` entry containing fields
│ ("nothing", "rpar", "Reflected Photosynthetically Active Radiation", "reflected_photosynthetically_active_radiation", "mol photons m^-2 s^-1", "The amount of photosynthetically active radiation reflected by leaves.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `tpar` entry containing fields
│ ("nothing", "tpar", "Transmitted Photosynthetically Active Radiation", "transmitted_photosynthetically_active_radiation", "mol photons m^-2 s^-1", "The amount of photosynthetically active radiation transmitted by leaves.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `lwn` entry containing fields
│ ("nothing", "lwn", "Net Longwave Radiation", "net_longwave_radiation", "W m^-2", "The net (down minus up) longwave radiation at the surface.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `swn` entry containing fields
│ ("nothing", "swn", "Net Shortwave Radiation", "net_shortwave_radiation", "W m^-2", "The net (down minus up) shortwave radiation at the surface.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `soc` entry containing fields
│ ("nothing", "soc", "Soil organic carbon", "soil_organic_carbon", "kg C m^-3", "Mass of organic carbon per volume of soil. (depth resolved)")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `airp` entry containing fields
│ ("nothing", "airp", "Air pressure", "air_pressure", "Pa", "The air pressure.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `rain` entry containing fields
│ ("nothing", "rain", "Rainfall", "rainfall", "m s^-1", "Precipitation of liquid water volume (m^3 of water per m^2 of ground per second).")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `lwd` entry containing fields
│ ("nothing", "lwd", "Down Longwave Radiation", "down_longwave_radiation", "W m^-2", "The downwelling longwave radiation at the surface.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `swd` entry containing fields
│ ("nothing", "swd", "Shortwave Radiation Downwards", "down_shortwave_radiation", "W m^-2", "The downwelling shortwave radiation at the surface.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `snow` entry containing fields
│ ("nothing", "snow", "Snowfall", "snowfall", "m s^-1", "The precipitation of snow in liquid water volume (m^3 of water per m^2 of ground per second).")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `precip` entry containing fields
│ ("nothing", "precip", "Total precipitation", "total_precipitation", "kg m^-2 s^-1", "The total flux from precipitation in kg of water per m^2 of ground per second).")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `tair` entry containing fields
│ ("nothing", "tair", "Air Temperature (K)", "tair", "K", "The air temperature at the lowest level of the atmosphere (coupled) or 2m level (prescribed).")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `qsfc` entry containing fields
│ ("nothing", "qsfc", "Surface Specific Humidity", "surface_specific_humidity", "", "Ratio of water vapor mass to total moist air parcel mass.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `ws` entry containing fields
│ ("nothing", "ws", "Wind Speed", "wind_speed", "m s^-1", "The average wind speed.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `infil` entry containing fields
│ ("nothing", "infil", "Infiltration", "infiltration", "m s^-1", "The flux of liquid water volume into the soil (m^3 of water per m^2 of ground per second).")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `salb` entry containing fields
│ ("nothing", "salb", "Soil Albedo", "surface albedo", "", "The mean of PAR and NIR albedo, which are calculated as α_soil_band = α_band_dry * (1 - S_e) + α_band_wet * S_e.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `shc` entry containing fields
│ ("nothing", "shc", "Soil Hydraulic Conductivity", "soil_hydraulic_conductivity", "m s^-1", "Soil hydraulic conductivity. (depth resolved)")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `stc` entry containing fields
│ ("nothing", "stc", "Soil Thermal Conductivity", "soil_thermal_conductivity", "W m^-1 K^-1", "Soil thermal conductivity. (depth resolved)")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `swp` entry containing fields
│ ("nothing", "swp", "Soil Water Potential", "soil_water_potential", "Pa", "Soil water potential. (depth resolved)")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `soilrn` entry containing fields
│ ("nothing", "soilrn", "Soil Net Radiation", "soil_net_radiation", "W m^-2", "Net radiation at the soil surface.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `soillhf` entry containing fields
│ ("nothing", "soillhf", "Soil Latent Heat Flux", "soil_Latent_Heat_Flux", "W m^-2", "Soil latent heat flux, the amount of liquid water evaporated by the soil, expressed in energy units (W m^-2).")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `soilshf` entry containing fields
│ ("nothing", "soilshf", "Soil Sensible Heat Flux", "soil_sensible_Heat_Flux", "W m^-2", "Soil sensible heat flux, the amount of energy exchanged between the soil and atmosphere to change the temperature of the soil.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `hr` entry containing fields
│ ("nothing", "hr", "Heterotrophic Respiration", "heterotrophic_respiration", "mol m^-2 s^-1", "The CO2 efflux at the soil surface due to microbial decomposition of soil organic matter. This is not necessarily equal to CO2 production by microbes, as co2 diffusion through the soil pores takes time.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `scd` entry containing fields
│ ("nothing", "scd", "Soil CO2 Diffusivity", "soil_co2_diffusivity", "m^2 s^-1", "The diffusivity of CO2 in the porous phase of the soil. Depends on soil texture, moisture, and temperature. (depth resolved)")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `scms` entry containing fields
│ ("nothing", "scms", "Soil CO2 Microbial Source", "soil_co2_microbial_source", "kg C m^-3 s^-1", "The production of CO2 by microbes in the soil. Vary by layers of soil depth. (depth resolved)")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `lwu` entry containing fields
│ ("nothing", "lwu", "Longwave Radiation Up", "longwave_radiation_up", "W m^-2", "Upwelling longwave radiation.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `swu` entry containing fields
│ ("nothing", "swu", "Shortwave Radiation Up", "shortwave_radiation_up", "W m^-2", "Upwelling shortwave radiation")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `et` entry containing fields
│ ("nothing", "et", "Evapotranspiration", "evapotranspiration", "kg m^-2 s^-1", "Total flux of water mass out of the surface.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `er` entry containing fields
│ ("nothing", "er", "Ecosystem Respiration", "ecosystem respiration", "mol CO2 m^-2 s^-1", "Total respiration flux out of the surface.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `sr` entry containing fields
│ ("nothing", "sr", "Surface Runoff", "surface_runoff", "m s^-1", "Water runoff at the surface, this is the water flowing horizontally above the ground.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `ssr` entry containing fields
│ ("nothing", "ssr", "Subsurface Runoff", "subsurface_runoff", "m s^-1", "Water runoff from below the surface")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `ct` entry containing fields
│ ("nothing", "ct", "Canopy Temperature", "canopy_temperature", "K", "Canopy temperature.")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `sco2` entry containing fields
│ ("nothing", "sco2", "Soil CO2", "soil_co2", "kg C m^3", "Concentration of CO2 in the porous air space of the soil. (depth resolved)")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `swc` entry containing fields
│ ("nothing", "swc", "Soil Water Content", "soil_water_content", "m^3 m^-3", "The volume of soil water per volume of soil. (depth resolved)")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `iwc` entry containing fields
│ ("nothing", "iwc", "Integrated Soil Water Mass in first 10cm", "soil_10cm_water_mass", "kg/m^2", "The integrated water mass to a depth of 10cm")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `si` entry containing fields
│ ("nothing", "si", "Soil Ice", "soil_ice", "m^3 m^-3", "The volume of soil ice per volume of soil. (depth resolved)")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `sie` entry containing fields
│ ("nothing", "sie", "Soil Internal Energy", "soil_internal_energy", "W m^-2", "The energy per volume of soil. (depth resolved)")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `swe` entry containing fields
│ ("nothing", "swe", "Snow water equivalent", "snow_water_equivalent", "m", "The height of liquid water if all snow melted")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `snd` entry containing fields
│ ("nothing", "snd", "Snow depth", "snow_depth", "m", "The snow depth")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Warning: overwriting diagnostic `snowc` entry containing fields
│ ("nothing", "snowc", "Snow cover fraction", "snow_cover_fraction", "", "The snow cover fraction")
└ @ ClimaLand.Diagnostics ~/work/ClimaLand.jl/ClimaLand.jl/src/diagnostics/diagnostic.jl:60
┌ Info: Progress
│   simulation_time = "1 week, 3 days"
│   n_steps_completed = 1000
│   wall_time_per_step = "0 seconds"
│   wall_time_total = "0 seconds"
│   wall_time_remaining = "0 seconds"
│   wall_time_spent = "0 seconds"
│   percent_complete = "2.8%"
│   estimated_sypd = "Inf"
│   date_now = 2025-08-13T20:43:57.776
└   estimated_finish_date = 2025-08-13T20:43:57.775
┌ Info: Progress
│   simulation_time = "2 weeks, 6 days"
│   n_steps_completed = 2000
│   wall_time_per_step = "0 seconds"
│   wall_time_total = "0 seconds"
│   wall_time_remaining = "0 seconds"
│   wall_time_spent = "0 seconds"
│   percent_complete = "5.7%"
│   estimated_sypd = "Inf"
│   date_now = 2025-08-13T20:44:05.845
└   estimated_finish_date = 2025-08-13T20:44:05.845
[ Info: Checking NaNs in bucket
[ Info: Checking NaNs in W
[ Info: Checking NaNs in T
[ Info: Checking NaNs in Ws
[ Info: Checking NaNs in σS
┌ Info: Progress
│   simulation_time = "4 weeks, 3 days"
│   n_steps_completed = 3000
│   wall_time_per_step = "15 milliseconds, 636 microseconds"
│   wall_time_total = "9 minutes, 9 seconds"
│   wall_time_remaining = "8 minutes, 22 seconds"
│   wall_time_spent = "46 seconds, 908 milliseconds"
│   percent_complete = "8.5%"
│   estimated_sypd = "157.586"
│   date_now = 2025-08-13T20:44:29.496
└   estimated_finish_date = 2025-08-13T20:52:52.495
┌ Info: Progress
│   simulation_time = "5 weeks, 6 days"
│   n_steps_completed = 4000
│   wall_time_per_step = "13 milliseconds, 822 microseconds"
│   wall_time_total = "8 minutes, 5 seconds"
│   wall_time_remaining = "7 minutes, 10 seconds"
│   wall_time_spent = "55 seconds, 289 milliseconds"
│   percent_complete = "11.4%"
│   estimated_sypd = "178.266"
│   date_now = 2025-08-13T20:44:37.680
└   estimated_finish_date = 2025-08-13T20:51:48.680
┌ Info: Progress
│   simulation_time = "7 weeks, 3 days"
│   n_steps_completed = 5000
│   wall_time_per_step = "12 milliseconds, 681 microseconds"
│   wall_time_total = "7 minutes, 25 seconds"
│   wall_time_remaining = "6 minutes, 22 seconds"
│   wall_time_spent = "1 minute, 3 seconds"
│   percent_complete = "14.2%"
│   estimated_sypd = "194.305"
│   date_now = 2025-08-13T20:44:45.797
└   estimated_finish_date = 2025-08-13T20:51:08.797
[ Info: Checking NaNs in bucket
[ Info: Checking NaNs in W
[ Info: Checking NaNs in T
[ Info: Checking NaNs in Ws
[ Info: Checking NaNs in σS
┌ Info: Progress
│   simulation_time = "8 weeks, 6 days"
│   n_steps_completed = 6000
│   wall_time_per_step = "11 milliseconds, 974 microseconds"
│   wall_time_total = "7 minutes, 753 milliseconds"
│   wall_time_remaining = "5 minutes, 48 seconds"
│   wall_time_spent = "1 minute, 11 seconds"
│   percent_complete = "17.1%"
│   estimated_sypd = "205.768"
│   date_now = 2025-08-13T20:44:54.240
└   estimated_finish_date = 2025-08-13T20:50:43.240
┌ Info: Progress
│   simulation_time = "10 weeks, 2 days"
│   n_steps_completed = 7000
│   wall_time_per_step = "11 milliseconds, 434 microseconds"
│   wall_time_total = "6 minutes, 41 seconds"
│   wall_time_remaining = "5 minutes, 21 seconds"
│   wall_time_spent = "1 minute, 20 seconds"
│   percent_complete = "19.9%"
│   estimated_sypd = "215.491"
│   date_now = 2025-08-13T20:45:02.433
└   estimated_finish_date = 2025-08-13T20:50:24.433
┌ Info: Progress
│   simulation_time = "11 weeks, 6 days"
│   n_steps_completed = 8000
│   wall_time_per_step = "11 milliseconds, 48 microseconds"
│   wall_time_total = "6 minutes, 28 seconds"
│   wall_time_remaining = "4 minutes, 59 seconds"
│   wall_time_spent = "1 minute, 28 seconds"
│   percent_complete = "22.8%"
│   estimated_sypd = "223.026"
│   date_now = 2025-08-13T20:45:10.777
└   estimated_finish_date = 2025-08-13T20:50:10.777
[ Info: Checking NaNs in bucket
[ Info: Checking NaNs in W
[ Info: Checking NaNs in T
[ Info: Checking NaNs in Ws
[ Info: Checking NaNs in σS
┌ Info: Progress
│   simulation_time = "13 weeks, 2 days"
│   n_steps_completed = 9000
│   wall_time_per_step = "10 milliseconds, 745 microseconds"
│   wall_time_total = "6 minutes, 17 seconds"
│   wall_time_remaining = "4 minutes, 40 seconds"
│   wall_time_spent = "1 minute, 36 seconds"
│   percent_complete = "25.6%"
│   estimated_sypd = "229.307"
│   date_now = 2025-08-13T20:45:19.102
└   estimated_finish_date = 2025-08-13T20:50:00.102
┌ Info: Progress
│   simulation_time = "14 weeks, 6 days"
│   n_steps_completed = 10000
│   wall_time_per_step = "10 milliseconds, 520 microseconds"
│   wall_time_total = "6 minutes, 9 seconds"
│   wall_time_remaining = "4 minutes, 24 seconds"
│   wall_time_spent = "1 minute, 45 seconds"
│   percent_complete = "28.5%"
│   estimated_sypd = "234.209"
│   date_now = 2025-08-13T20:45:27.599
└   estimated_finish_date = 2025-08-13T20:49:52.599
┌ Info: Progress
│   simulation_time = "16 weeks, 2 days"
│   n_steps_completed = 11000
│   wall_time_per_step = "10 milliseconds, 342 microseconds"
│   wall_time_total = "6 minutes, 3 seconds"
│   wall_time_remaining = "4 minutes, 9 seconds"
│   wall_time_spent = "1 minute, 53 seconds"
│   percent_complete = "31.3%"
│   estimated_sypd = "238.257"
│   date_now = 2025-08-13T20:45:36.153
└   estimated_finish_date = 2025-08-13T20:49:46.153
[ Info: Checking NaNs in bucket
[ Info: Checking NaNs in W
[ Info: Checking NaNs in T
[ Info: Checking NaNs in Ws
[ Info: Checking NaNs in σS
┌ Info: Progress
│   simulation_time = "17 weeks, 6 days"
│   n_steps_completed = 12000
│   wall_time_per_step = "10 milliseconds, 198 microseconds"
│   wall_time_total = "5 minutes, 58 seconds"
│   wall_time_remaining = "3 minutes, 55 seconds"
│   wall_time_spent = "2 minutes, 2 seconds"
│   percent_complete = "34.2%"
│   estimated_sypd = "241.618"
│   date_now = 2025-08-13T20:45:44.769
└   estimated_finish_date = 2025-08-13T20:49:40.769
┌ Info: Progress
│   simulation_time = "19 weeks, 2 days"
│   n_steps_completed = 13000
│   wall_time_per_step = "10 milliseconds, 82 microseconds"
│   wall_time_total = "5 minutes, 54 seconds"
│   wall_time_remaining = "3 minutes, 43 seconds"
│   wall_time_spent = "2 minutes, 11 seconds"
│   percent_complete = "37.0%"
│   estimated_sypd = "244.391"
│   date_now = 2025-08-13T20:45:53.463
└   estimated_finish_date = 2025-08-13T20:49:37.462
┌ Info: Progress
│   simulation_time = "20 weeks, 5 days"
│   n_steps_completed = 14000
│   wall_time_per_step = "9 milliseconds, 987 microseconds"
│   wall_time_total = "5 minutes, 50 seconds"
│   wall_time_remaining = "3 minutes, 31 seconds"
│   wall_time_spent = "2 minutes, 19 seconds"
│   percent_complete = "39.8%"
│   estimated_sypd = "246.704"
│   date_now = 2025-08-13T20:46:02.221
└   estimated_finish_date = 2025-08-13T20:49:34.221
[ Info: Checking NaNs in bucket
[ Info: Checking NaNs in W
[ Info: Checking NaNs in T
[ Info: Checking NaNs in Ws
[ Info: Checking NaNs in σS
┌ Info: Progress
│   simulation_time = "22 weeks, 2 days"
│   n_steps_completed = 15000
│   wall_time_per_step = "9 milliseconds, 910 microseconds"
│   wall_time_total = "5 minutes, 48 seconds"
│   wall_time_remaining = "3 minutes, 19 seconds"
│   wall_time_spent = "2 minutes, 28 seconds"
│   percent_complete = "42.7%"
│   estimated_sypd = "248.641"
│   date_now = 2025-08-13T20:46:11.043
└   estimated_finish_date = 2025-08-13T20:49:31.043
┌ Info: Progress
│   simulation_time = "23 weeks, 5 days"
│   n_steps_completed = 16000
│   wall_time_per_step = "9 milliseconds, 841 microseconds"
│   wall_time_total = "5 minutes, 45 seconds"
│   wall_time_remaining = "3 minutes, 8 seconds"
│   wall_time_spent = "2 minutes, 37 seconds"
│   percent_complete = "45.5%"
│   estimated_sypd = "250.371"
│   date_now = 2025-08-13T20:46:19.857
└   estimated_finish_date = 2025-08-13T20:49:28.857
┌ Info: Progress
│   simulation_time = "25 weeks, 2 days"
│   n_steps_completed = 17000
│   wall_time_per_step = "9 milliseconds, 783 microseconds"
│   wall_time_total = "5 minutes, 43 seconds"
│   wall_time_remaining = "2 minutes, 57 seconds"
│   wall_time_spent = "2 minutes, 46 seconds"
│   percent_complete = "48.4%"
│   estimated_sypd = "251.866"
│   date_now = 2025-08-13T20:46:28.705
└   estimated_finish_date = 2025-08-13T20:49:26.705
[ Info: Checking NaNs in bucket
[ Info: Checking NaNs in W
[ Info: Checking NaNs in T
[ Info: Checking NaNs in Ws
[ Info: Checking NaNs in σS
┌ Info: Progress
│   simulation_time = "26 weeks, 5 days"
│   n_steps_completed = 18000
│   wall_time_per_step = "9 milliseconds, 733 microseconds"
│   wall_time_total = "5 minutes, 41 seconds"
│   wall_time_remaining = "2 minutes, 46 seconds"
│   wall_time_spent = "2 minutes, 55 seconds"
│   percent_complete = "51.2%"
│   estimated_sypd = "253.154"
│   date_now = 2025-08-13T20:46:37.593
└   estimated_finish_date = 2025-08-13T20:49:24.593
┌ Info: Progress
│   simulation_time = "28 weeks, 1 day"
│   n_steps_completed = 19000
│   wall_time_per_step = "9 milliseconds, 690 microseconds"
│   wall_time_total = "5 minutes, 40 seconds"
│   wall_time_remaining = "2 minutes, 36 seconds"
│   wall_time_spent = "3 minutes, 4 seconds"
│   percent_complete = "54.1%"
│   estimated_sypd = "254.271"
│   date_now = 2025-08-13T20:46:46.514
└   estimated_finish_date = 2025-08-13T20:49:23.514
┌ Info: Progress
│   simulation_time = "29 weeks, 5 days"
│   n_steps_completed = 20000
│   wall_time_per_step = "9 milliseconds, 650 microseconds"
│   wall_time_total = "5 minutes, 39 seconds"
│   wall_time_remaining = "2 minutes, 26 seconds"
│   wall_time_spent = "3 minutes, 13 seconds"
│   percent_complete = "56.9%"
│   estimated_sypd = "255.319"
│   date_now = 2025-08-13T20:46:55.409
└   estimated_finish_date = 2025-08-13T20:49:22.409
[ Info: Checking NaNs in bucket
[ Info: Checking NaNs in W
[ Info: Checking NaNs in T
[ Info: Checking NaNs in Ws
[ Info: Checking NaNs in σS
┌ Info: Progress
│   simulation_time = "31 weeks, 1 day"
│   n_steps_completed = 21000
│   wall_time_per_step = "9 milliseconds, 615 microseconds"
│   wall_time_total = "5 minutes, 37 seconds"
│   wall_time_remaining = "2 minutes, 15 seconds"
│   wall_time_spent = "3 minutes, 21 seconds"
│   percent_complete = "59.8%"
│   estimated_sypd = "256.252"
│   date_now = 2025-08-13T20:47:04.322
└   estimated_finish_date = 2025-08-13T20:49:20.322
┌ Info: Progress
│   simulation_time = "32 weeks, 5 days"
│   n_steps_completed = 22000
│   wall_time_per_step = "9 milliseconds, 580 microseconds"
│   wall_time_total = "5 minutes, 36 seconds"
│   wall_time_remaining = "2 minutes, 5 seconds"
│   wall_time_spent = "3 minutes, 30 seconds"
│   percent_complete = "62.6%"
│   estimated_sypd = "257.183"
│   date_now = 2025-08-13T20:47:13.172
└   estimated_finish_date = 2025-08-13T20:49:19.172
┌ Info: Progress
│   simulation_time = "34 weeks, 1 day"
│   n_steps_completed = 23000
│   wall_time_per_step = "9 milliseconds, 546 microseconds"
│   wall_time_total = "5 minutes, 35 seconds"
│   wall_time_remaining = "1 minute, 55 seconds"
│   wall_time_spent = "3 minutes, 39 seconds"
│   percent_complete = "65.5%"
│   estimated_sypd = "258.115"
│   date_now = 2025-08-13T20:47:21.957
└   estimated_finish_date = 2025-08-13T20:49:17.957
[ Info: Checking NaNs in bucket
[ Info: Checking NaNs in W
[ Info: Checking NaNs in T
[ Info: Checking NaNs in Ws
[ Info: Checking NaNs in σS
┌ Info: Progress
│   simulation_time = "35 weeks, 5 days"
│   n_steps_completed = 24000
│   wall_time_per_step = "9 milliseconds, 511 microseconds"
│   wall_time_total = "5 minutes, 34 seconds"
│   wall_time_remaining = "1 minute, 45 seconds"
│   wall_time_spent = "3 minutes, 48 seconds"
│   percent_complete = "68.3%"
│   estimated_sypd = "259.052"
│   date_now = 2025-08-13T20:47:30.675
└   estimated_finish_date = 2025-08-13T20:49:16.675
┌ Info: Progress
│   simulation_time = "37 weeks, 1 day"
│   n_steps_completed = 25000
│   wall_time_per_step = "9 milliseconds, 476 microseconds"
│   wall_time_total = "5 minutes, 32 seconds"
│   wall_time_remaining = "1 minute, 36 seconds"
│   wall_time_spent = "3 minutes, 56 seconds"
│   percent_complete = "71.2%"
│   estimated_sypd = "260.032"
│   date_now = 2025-08-13T20:47:39.290
└   estimated_finish_date = 2025-08-13T20:49:16.290
┌ Info: Progress
│   simulation_time = "38 weeks, 4 days"
│   n_steps_completed = 26000
│   wall_time_per_step = "9 milliseconds, 439 microseconds"
│   wall_time_total = "5 minutes, 31 seconds"
│   wall_time_remaining = "1 minute, 26 seconds"
│   wall_time_spent = "4 minutes, 5 seconds"
│   percent_complete = "74.0%"
│   estimated_sypd = "261.05"
│   date_now = 2025-08-13T20:47:47.806
└   estimated_finish_date = 2025-08-13T20:49:14.806
[ Info: Checking NaNs in bucket
[ Info: Checking NaNs in W
[ Info: Checking NaNs in T
[ Info: Checking NaNs in Ws
[ Info: Checking NaNs in σS
┌ Info: Progress
│   simulation_time = "40 weeks, 1 day"
│   n_steps_completed = 27000
│   wall_time_per_step = "9 milliseconds, 395 microseconds"
│   wall_time_total = "5 minutes, 30 seconds"
│   wall_time_remaining = "1 minute, 16 seconds"
│   wall_time_spent = "4 minutes, 13 seconds"
│   percent_complete = "76.8%"
│   estimated_sypd = "262.256"
│   date_now = 2025-08-13T20:47:56.073
└   estimated_finish_date = 2025-08-13T20:49:13.073
┌ Info: Progress
│   simulation_time = "41 weeks, 4 days"
│   n_steps_completed = 28000
│   wall_time_per_step = "9 milliseconds, 353 microseconds"
│   wall_time_total = "5 minutes, 28 seconds"
│   wall_time_remaining = "1 minute, 6 seconds"
│   wall_time_spent = "4 minutes, 21 seconds"
│   percent_complete = "79.7%"
│   estimated_sypd = "263.431"
│   date_now = 2025-08-13T20:48:04.295
└   estimated_finish_date = 2025-08-13T20:49:11.295
┌ Info: Progress
│   simulation_time = "43 weeks, 1 day"
│   n_steps_completed = 29000
│   wall_time_per_step = "9 milliseconds, 308 microseconds"
│   wall_time_total = "5 minutes, 27 seconds"
│   wall_time_remaining = "57 seconds, 116 milliseconds"
│   wall_time_spent = "4 minutes, 29 seconds"
│   percent_complete = "82.5%"
│   estimated_sypd = "264.711"
│   date_now = 2025-08-13T20:48:12.337
└   estimated_finish_date = 2025-08-13T20:49:10.337
[ Info: Checking NaNs in bucket
[ Info: Checking NaNs in W
[ Info: Checking NaNs in T
[ Info: Checking NaNs in Ws
[ Info: Checking NaNs in σS
┌ Info: Progress
│   simulation_time = "44 weeks, 4 days"
│   n_steps_completed = 30000
│   wall_time_per_step = "9 milliseconds, 268 microseconds"
│   wall_time_total = "5 minutes, 25 seconds"
│   wall_time_remaining = "47 seconds, 601 milliseconds"
│   wall_time_spent = "4 minutes, 38 seconds"
│   percent_complete = "85.4%"
│   estimated_sypd = "265.862"
│   date_now = 2025-08-13T20:48:20.437
└   estimated_finish_date = 2025-08-13T20:49:08.437
┌ Info: Progress
│   simulation_time = "46 weeks, 22 hours"
│   n_steps_completed = 31000
│   wall_time_per_step = "9 milliseconds, 227 microseconds"
│   wall_time_total = "5 minutes, 24 seconds"
│   wall_time_remaining = "38 seconds, 166 milliseconds"
│   wall_time_spent = "4 minutes, 46 seconds"
│   percent_complete = "88.2%"
│   estimated_sypd = "267.023"
│   date_now = 2025-08-13T20:48:28.456
└   estimated_finish_date = 2025-08-13T20:49:07.456
┌ Info: Progress
│   simulation_time = "47 weeks, 4 days"
│   n_steps_completed = 32000
│   wall_time_per_step = "9 milliseconds, 184 microseconds"
│   wall_time_total = "5 minutes, 22 seconds"
│   wall_time_remaining = "28 seconds, 802 milliseconds"
│   wall_time_spent = "4 minutes, 53 seconds"
│   percent_complete = "91.1%"
│   estimated_sypd = "268.283"
│   date_now = 2025-08-13T20:48:36.297
└   estimated_finish_date = 2025-08-13T20:49:05.297
[ Info: Checking NaNs in bucket
[ Info: Checking NaNs in W
[ Info: Checking NaNs in T
[ Info: Checking NaNs in Ws
[ Info: Checking NaNs in σS
┌ Info: Progress
│   simulation_time = "49 weeks, 18 hours"
│   n_steps_completed = 33000
│   wall_time_per_step = "9 milliseconds, 148 microseconds"
│   wall_time_total = "5 minutes, 21 seconds"
│   wall_time_remaining = "19 seconds, 540 milliseconds"
│   wall_time_spent = "5 minutes, 1 second"
│   percent_complete = "93.9%"
│   estimated_sypd = "269.347"
│   date_now = 2025-08-13T20:48:44.284
└   estimated_finish_date = 2025-08-13T20:49:04.284
┌ Info: Progress
│   simulation_time = "50 weeks, 4 days"
│   n_steps_completed = 34000
│   wall_time_per_step = "9 milliseconds, 113 microseconds"
│   wall_time_total = "5 minutes, 20 seconds"
│   wall_time_remaining = "10 seconds, 352 milliseconds"
│   wall_time_spent = "5 minutes, 9 seconds"
│   percent_complete = "96.8%"
│   estimated_sypd = "270.384"
│   date_now = 2025-08-13T20:48:52.239
└   estimated_finish_date = 2025-08-13T20:49:03.239
┌ Info: Progress
│   simulation_time = "52 weeks, 14 hours"
│   n_steps_completed = 35000
│   wall_time_per_step = "9 milliseconds, 76 microseconds"
│   wall_time_total = "5 minutes, 18 seconds"
│   wall_time_remaining = "1 second, 234 milliseconds"
│   wall_time_spent = "5 minutes, 17 seconds"
│   percent_complete = "99.6%"
│   estimated_sypd = "271.48"
│   date_now = 2025-08-13T20:49:00.065
└   estimated_finish_date = 2025-08-13T20:49:02.065
[ Info: Checking NaNs in bucket
[ Info: Checking NaNs in W
[ Info: Checking NaNs in T
[ Info: Checking NaNs in Ws
[ Info: Checking NaNs in σS

Make some plots:

short_names = ["tsfc", "lhf", "shf", "wsoil"]

LandSimVis.make_annual_timeseries(
    simulation;
    savedir = ".",
    short_names,
    plot_name = "bucket_annual_timeseries.pdf",
)

LandSimVis.make_heatmaps(
    simulation;
    savedir = ".",
    short_names,
    date = stop_date,
    plot_name = "bucket_figures.pdf",
)


This page was generated using Literate.jl.