Drivers

ClimaLand.PrescribedAtmosphereType
PrescribedAtmosphere{FT, CA, DT} <: AbstractAtmosphericDrivers{FT}

Container for holding prescribed atmospheric drivers and other information needed for computing turbulent surface fluxes when driving land models in standalone mode.

The default CO2 concentration is a constant as a function of time, equal to 4.2e-4 mol/mol.

Since not all models require co2 concentration, the default for that is nothing.

  • liquid_precip: Precipitation (m/s) function of time: positive by definition

  • snow_precip: Snow precipitation (m/s) function of time: positive by definition

  • T: Prescribed atmospheric temperature (function of time) at the reference height (K)

  • u: Prescribed wind speed (function of time) at the reference height (m/s)

  • q: Prescribed specific humidity (function of time) at the reference height (_)

  • P: Prescribed air pressure (function of time) at the reference height (Pa)

  • c_co2: CO2 concentration in atmosphere (mol/mol)

  • start_date: Start date - the datetime corresponding to t=0 for the simulation

  • h: Reference height (m), relative to surface elevation

  • gustiness: Minimum wind speed (gustiness; m/s)

  • thermo_params: Thermodynamic parameters

source
ClimaLand.PrescribedPrecipitationType
PrescribedPrecipitation{FT, LP} <: AbstractAtmosphericDrivers{FT}

Container for holding prescribed precipitation driver for models which only require precipitation (RichardsModel).

  • liquid_precip: Precipitation (m/s) function of time: positive by definition
source
ClimaLand.PrescribedRadiativeFluxesType
PrescribedRadiativeFluxes{FT, SW, LW, DT, T, TP} <: AbstractRadiativeDrivers{FT}

Container for the prescribed radiation functions needed to drive land models in standalone mode.

Note that some models require the zenith angle AND diffuse fraction. The diffuse fraction may be provided directly (of type TimeVaryingInput), or it may be computed empirically. In the latter case, it requires the thermodynamic parameters as well to compute.

Therefore, the allowed combinations are:

  1. Zenith angle and diffuse fraction not needed: zenith angle=nothing, thermo_params=nothing, diffuse fraction=nothing
  2. Zenith angle provided and diffuse fraction computed empirically: thermo params used, diffuse fraction=nothing
  3. Zenith angle provided and diffuse fraction provided: thermo_params not used, diffuse fraction TimeVaryingInput
  • SW_d: Downward shortwave radiation function of time (W/m^2): positive indicates towards surface

  • frac_diff: Diffuse Fraction of shortwave radiation (unitless, [0,1])

  • LW_d: Downward longwave radiation function of time (W/m^2): positive indicates towards surface

  • start_date: Start date - the datetime corresponding to t=0 for the simulation

  • θs: Sun zenith angle, in radians

  • thermo_params: Thermodynamic parameters

source
ClimaLand.CoupledAtmosphereType
CoupledAtmosphere{FT} <: AbstractAtmosphericDrivers{FT}

To be used when coupling to an atmosphere model. Contains fields that are used to compute surface fluxes in the coupled setup.

When constructed without a space, the struct doesn't contain anything, but it still acts as a flag that fluxes have been updated by the coupler and don't need to be recomputed. When constructed with a space, the struct contains the fields needed to compute surface fluxes in the coupled setup, which are accessed by ClimaCoupler.

source
ClimaLand.CoupledRadiativeFluxesType
CoupledRadiativeFluxes{
    FT,
    F <: Union{Function, Nothing},
    T,
} <: AbstractRadiativeDrivers{FT}

To be used when coupling to an atmosphere model. Either both θs and start_date must be nothing, or both must not be nothing`.

During the driver update, cosθs is unchanged if θs is nothing. This behavior differs from the PrescribedRadiativeFluxes where cosθs set to NaN if θs is nothing. Otherwise, θs recieves the following arguments: (time_from_start, start_date), and is expected to return zenith angle at the given time.

  • θs: Function that fills a climacore field with the zenith angle given the following arguments: (timefromstart, start_date)

  • start_date: Start date - the datetime corresponding to t=0 for the simulation

source
ClimaLand.PrescribedGroundConditionsType
 PrescribedGroundConditions <: AbstractGroundConditions

A container for holding prescribed ground conditions needed by the canopy model when running the canopy in standalone mode, including the surface temperature, albedo, and emissivity, and soil water content, porosity, residual water fraction, and hydrology closure model.

Note that internally we enforce that the soil water content θ must be within the range (θ\_r, ν].

  • θ: Prescribed soil water content (m^3/m^3) in the root zone as a function of time

  • T: Prescribed ground surface temperature (K) as a function of time

  • α_PAR: Ground albedo for PAR

  • α_NIR: Ground albedo for NIR

  • ϵ: Ground emissivity

  • ν: Soil porosity (m^3/m^3

  • θ_r: The soil residual water fraction (m^3/m^3)

  • hydrology_cm: The soil hydrology closure model: van Genuchten or Brooks and Corey

source
ClimaLand.PrognosticGroundConditionsType
 PrognosticGroundConditions <: Canopy.AbstractGroundConditions

A type of AbstractGroundConditions to use when running the CanopyModel as part of an integrated model, i.e. with a prognostic EnergyHydrology soil model and either with or without a snow model.

Note that this struct is linked with the EnergyHydrology/SnowModel models. If we ever had a different soil model, we might need to construct a different PrognosticGroundConditions because the fields may be stored in different places.

source
ClimaLand.turbulent_fluxes!Function
turbulent_fluxes!(dest,
                  atmos::PrescribedAtmosphere,
                  model::AbstractModel,
                  Y::ClimaCore.Fields.FieldVector,
                  p::NamedTuple,
                  t
                  )

Computes the turbulent surface flux terms at the ground for a standalone simulation, including turbulent energy fluxes as well as the water vapor flux (in units of m^3/m^2/s of water). Positive fluxes indicate flow from the ground to the atmosphere.

It solves for these given atmospheric conditions, stored in atmos, model parameters, and the surface conditions.

source
turbulent_fluxes!(dest,
                atmos::CoupledAtmosphere,
                model::AbstractModel,
                Y,
                p,
                t)

Computes the turbulent surface fluxes terms at the ground for a coupled simulation. In this case, the coupler has already computed turbulent fluxes and updated them in each of the component models, so this function does nothing.

source
function ClimaLand.turbulent_fluxes!(
    dest,
    atmos::PrescribedAtmosphere,
    model::CanopyModel,
    Y::ClimaCore.Fields.FieldVector,
    p::NamedTuple,
    t,
)

A canopy specific function for compute turbulent fluxes with the atmosphere; returns the latent heat flux, sensible heat flux, vapor flux, and aerodynamic resistance.

We cannot use the default version in src/shared_utilities/drivers.jl because the canopy requires a different resistance for vapor and sensible heat fluxes, and the resistances depend on ustar, which we must compute using SurfaceFluxes before adjusting to account for these resistances.

source
ClimaLand.turbulent_fluxes_at_a_pointFunction
turbulent_fluxes_at_a_point(return_extra_fluxes, args...)

This is a wrapper function that allows us to dispatch on the type of return_extra_fluxes as we compute the turbulent fluxes pointwise. This is needed because space for the extra fluxes is only allocated in the cache when running with a CoupledAtmosphere. The function compute_turbulent_fluxes_at_a_point does the actual flux computation.

The return_extra_fluxes argument indicates whether to return the following:

  • momentum fluxes (ρτxz, ρτyz)
  • buoyancy flux (buoy_flux)
source
ClimaLand.set_atmos_ts!Function
set_atmos_ts!(ts_in, atmos::PrescribedAtmosphere{FT}, p)

Fill the pre-allocated ts_in Field with a thermodynamic state computed from the atmosphere.

source
ClimaLand.surface_air_densityFunction
ClimaLand.surface_air_density(
                atmos::CoupledAtmosphere,
                model::BucketModel,
                Y,
                p,
                _...,
            )

Returns the air density at the surface in the case of a coupled simulation.

This requires the field ρ_sfc to be present in the cache p under the name of the model.

source
surface_air_density(
                    atmos::PrescribedAtmosphere,
                    model::AbstractModel,
                    Y,
                    p,
                    t,
                    T_sfc,
                    )

A helper function which returns the surface air density; this assumes that the model has a property called parameters containing earth_param_set.

We additionally include the atmos type as an argument because the surface air density computation will change between a coupled simulation and a prescibed atmos simulation.

Extending this function for your model is only necessary if you need to compute the air density in a different way.

source
surface_air_density(
                    atmos::CoupledAtmosphere,
                    model::AbstractModel,
                    Y,
                    p,
                    t,
                    T_sfc,
                    )

A helper function which returns the surface air density; this assumes that the model has a property called parameters containing earth_param_set.

This method is similar to the general method above, except in this case we get the thermodynamic parameters from the atmos object. This is used when running with a coupled atmosphere.

source
ClimaLand.surface_temperatureFunction
ClimaLand.surface_temperature(
    model::EnergyHydrology{FT},
    Y,
    p,
    t,
) where {FT}

Returns the surface temperature field of the EnergyHydrology soil model.

The assumption is that the soil surface temperature is the same as the temperature at the center of the first soil layer.

source
ClimaLand.surface_temperature(model::SnowModel, Y, p)

a helper function which returns the surface temperature for the snow model, which is stored in the aux state.

source
ClimaLand.surface_temperature(model::BucketModel, Y, p)

a helper function which returns the surface temperature for the bucket model, which is stored in the aux state.

source
surface_temperature(model::AbstractModel, Y, p, t)

A helper function which returns the surface temperature for a given model, needed because different models compute and store surface temperature in different ways and places.

Extending this function for your model is only necessary if you need to compute surface fluxes and radiative fluxes at the surface using the functions in this file.

source
ClimaLand.surface_temperature(model::CanopyModel, Y, p, t)

A helper function which returns the temperature for the canopy model.

source
ClimaLand.surface_resistanceFunction
surface_resistance(model::AbstractModel, Y, p, t)

A helper function which returns the surface resistance for a given model, needed because different models compute and store surface resistance in different ways and places.

Extending this function for your model is only necessary if you need to compute surface fluxes and radiative fluxes at the surface using the functions in this file.

The default is 0, which is no additional resistance aside from the usual aerodynamic resistance from MOST.

source
ClimaLand.surface_resistance(
    model::CanopyModel{FT},
    Y,
    p,
    t,
) where {FT}

Returns the stomatal resistance field of the CanopyModel canopy.

source
ClimaLand.surface_specific_humidityFunction
ClimaLand.surface_specific_humidity(atmos::PrescribedAtmosphere, model::SnowModel, Y, p, _...)

Returns the precomputed specific humidity over snow as a weighted fraction of the saturated specific humidity over liquid and frozen water.

In the case of a PrescibedAtmoshere, the atmosphere thermal state is accessed from the drivers in the cache.

source
ClimaLand.surface_specific_humidity(atmos::CoupledAtmosphere, model::SnowModel, Y, p, _...)

Returns the precomputed specific humidity over snow as a weighted fraction of the saturated specific humidity over liquid and frozen water.

In the case of a CoupledAtmoshere, the atmosphere thermal state is accessed from the atmosphere object.

source
ClimaLand.surface_specific_humidity(atmos, model::BucketModel, Y, p)

a helper function which returns the surface specific humidity for the bucket model, which is stored in the aux state.

source
surface_specific_humidity(atmos, model::AbstractModel, Y, p, T_sfc, ρ_sfc)

A helper function which returns the surface specific humidity for a given model, needed because different models compute and store q_sfc in different ways and places.

Extending this function for your model is only necessary if you need to compute surface fluxes and radiative fluxes at the surface using the functions in this file.

source
ClimaLand.displacement_heightFunction
displacement_height(model::AbstractModel, Y, p)

A helper function which returns the displacement height for a given model; the default is zero.

Extending this function for your model is only necessary if you need to compute surface fluxes and radiative fluxes at the surface using the functions in this file.

source
ClimaLand.displacment_height(model::CanopyModel, Y, p)

A helper function which returns the displacement height for the canopy model.

See Cowan 1968; Brutsaert 1982, pp. 113–116; Campbell and Norman 1998, p. 71; Shuttleworth 2012, p. 343; Monteith and Unsworth 2013, p. 304.

source
ClimaLand.default_zenith_angleFunction
default_zenith_angle(
    t::T,
    start_date::Dates.DateTime;
    latitude::LT,
    longitude::LT,
    insol_params::Insolation.Parameters.InsolationParameters{FT},
)

Calculate zenith angle with Insolation for the given start date, insolation parameters, latitude, and longitude.

latitude and longitude can be a collections or a Number.

source
ClimaLand.prescribed_forcing_era5Function
prescribed_forcing_era5(start_date,
                        stop_date,
                        surface_space,
                        toml_dict::CP.ParamDict,
                        FT;
                        use_lowres_forcing = false,
                        gustiness=1,
                        max_wind_speed = nothing,
                        c_co2 = TimeVaryingInput((t) -> 4.2e-4),
                        time_interpolation_method = LinearInterpolation(PeriodicCalendar()),
                        regridder_type = :InterpolationsRegridder,
                        context = nothing,
                        )

A helper function which constructs the PrescribedAtmosphere and PrescribedRadiativeFluxes from ERA5 data stored in NetCDF files.

There are two versions of the ERA5 data available through ClimaLand:

  • a high resolution version (1° x 1°) available for the years 1979-2024
  • a low resolution version (8° x 8°) available only for the year 2008

The high resolution version will be used if use_lowres_forcing = false (the default). This artifact is used for global runs on compute clusters, but is too large to be used for local testing, and requires you to have acquired the data in advance. The low resolution version will be used if use_lowres_forcing = true. If the simulation dates are outside of 2008, the 2008 data will be reused for each year of simulation. This artifact is recommended for local testing or quick runs where accuracy is less critical.

The method for temporal interpolation is controlled via the time_interpolation_method kwarg. We suggest LinearInterpolation(PeriodicCalendar()), which implies linear interpolation in time; the inner argument implies how extrapolation outside the bounds of the data is handled. For example, the ERA5 forcing data we use is hourly, which implies Dec 31 of the last year of the data, at midnight, is not in the data. With the PeriodicCalendar() option, the interpolated value in the data at Jan 1 at timestamp 00 of the first year of the simulation will be used. For the low-resolution forcing data, which only exists for 2008, multi-year runs will repeat the forcing. If this behavior is not what you want, you can change the extrapolation argument to error LinearInterpolation(Throw()) or extrapolate by using the last value `LinearInterpolation(Flat()). More information is available in the ClimaUtilities documentation: https://clima.github.io/ClimaUtilities.jl/dev/inputs/#Extrapolation-boundary-conditions.

The ClimaLand default is to use nearest neighbor spatial interpolation for low resolution forcing, and linear spatial interpolation for high resolution forcing.

Clipped values

High wind speed anomalies (10-100x increase and decrease over a period of a several hours) appear in the ERA5 reanalysis data. These generate very large surface fluxes (due to wind speeds up to 300 m/s), which lead to instability. The kwarg max_wind_speed, with a value give in m/s, is used to clip these if it is not nothing. See here.

Full high resolution dataset available on clima cluster only

The full 40 year dataset of high resolution ERA5 data is only available on the clima cluster.

source
ClimaLand.prescribed_perturbed_temperature_era5Function
prescribed_perturbed_temperature_era5(era5_ncdata_path,
                         surface_space,
                         start_date,
                         toml_dict::CP.ParamDict,
                         ΔT,
                         FT;
                         gustiness=1,
                         max_wind_speed = nothing,
                         c_co2 = TimeVaryingInput((t) -> 4.2e-4),
                         time_interpolation_method = LinearInterpolation(PeriodicCalendar(Dates.Year(1), DateTime(Dates.year(stop_date)))),
                         regridder_type = :InterpolationsRegridder,
                         interpolation_method = Interpolations.Constant(),)

A helper function which constructs the PrescribedAtmosphere and PrescribedRadiativeFluxes from a file path pointing to the ERA5 data in a netcdf file, the surfacespace, the start date, and the `tomldict`, applying a change in the instantaneous temperature at each point of ΔT, while keeping the relative humidity fixed. The LWd shifts as LWd -> LW_d + aΔT, with a= 2W/m^2/K a surface climate sensitivity parameter.

Please see the documentation for prescribed_forcing_era5 for more information.

source
ClimaLand.prescribed_perturbed_rh_era5Function
 prescribed_perturbed_rh_era5(era5_ncdata_path,
                         surface_space,
                         start_date,
                         toml_dict::CP.ParamDict,
                         Δrh,
                         FT;
                         gustiness=1,
                         max_wind_speed = nothing,
                         c_co2 = TimeVaryingInput((t) -> 4.2e-4),
                         time_interpolation_method = LinearInterpolation(PeriodicCalendar(Dates.Year(1), DateTime(Dates.year(stop_date)))),
                         regridder_type = :InterpolationsRegridder,
                         interpolation_method = Interpolations.Constant(),)

A helper function which constructs the PrescribedAtmosphere and PrescribedRadiativeFluxes from a file path pointing to the ERA5 data in a netcdf file, the surfacespace, the start date, and the `tomldict`, applying a change in the instantaneous change to relative humidity at each point of Δrh. The perturbed rh is clipped to be within the range (0,1].

Please see the documentation for prescribed_forcing_era5 for more information.

source
ClimaLand.prescribed_analytic_forcingFunction
 prescribed_analytic_forcing(FT = Float32;
                             toml_dict::CP.ParamDict,
                             start_date = DateTime(2005),
                             SW_d = (t) -> 0,
                             LW_d = (t) -> 5.67e-8 * 280.0^4.0,
                             precip = (t) -> 0, # no precipitation
                             T_atmos = (t) -> 280.0,
                             u_atmos = (t) -> 1.0,
                             q_atmos = (t) -> 0.0, # no atmos water
                             h_atmos = FT(1e-8),
                             P_atmos = (t) -> 101325,
                             atmos = PrescribedAtmosphere(
                                 TimeVaryingInput(precip),
                                 TimeVaryingInput(precip),
                                 TimeVaryingInput(T_atmos),
                                 TimeVaryingInput(u_atmos),
                                 TimeVaryingInput(q_atmos),
                                 TimeVaryingInput(P_atmos),
                                 start_date,
                                 h_atmos,
                                 LP.LandParameters(toml_dict),
                             ),
                             radiation = PrescribedRadiativeFluxes(
                                 FT,
                                 TimeVaryingInput(SW_d),
                                 TimeVaryingInput(LW_d),
                                 start_date,
                             ),
                         )

A helper function which constructs the PrescribedAtmosphere and PrescribedRadiativeFluxes for a simple analytic case.

source
ClimaLand.net_radiation!Function
net_radiation!(dest::ClimaCore.Fields.Field,
               radiation::PrescribedRadiativeFluxes{FT},
               model::AbstractModel{FT},
               Y::ClimaCore.Fields.FieldVector,
               p::NamedTuple,
               t,
               ) where {FT}

Computes net radiative fluxes for a prescribed incoming longwave and shortwave radiation.

source
net_radiation!(dest,
              radiation::CoupledRadiativeFluxes,
              model::AbstractModel,
              Y,
              p,
              t)

Computes the net radiative flux at the ground for a coupled simulation. Your model cache must contain the field R_n.

source
ClimaLand.vapor_pressure_deficitFunction
vapor_pressure_deficit(T_air, P_air, q_air, thermo_params)

Computes the vapor pressure deficit for air with temperature Tair, pressure Pair, and specific humidity qair, using thermoparams, a Thermodynamics.jl param set.

source
ClimaLand.get_driversFunction
ClimaLand.get_drivers(model::RichardsModel)

Returns the driver variable symbols for the RichardsModel; these depend on the boundary condition type and currently only are required for the RichardsAtmosDrivenFluxBC, which is driven by a prescribed time and space varying precipitation.

source
ClimaLand.get_drivers(model::SnowModel)

Returns the driver variable symbols for the SnowModel.

source
get_drivers(model::AbstractModel)

Returns the driver objects for the model - atmospheric and radiative forcing, etc - as a tuple (atmos, radiation, ...). If no drivers are needed by a model, an empty tuple should be returned

source
ClimaLand.get_drivers(model::LandModel)

Returns the "drivers", or forcing variables, for the LandModel.

These consist of atmospheric and radiative forcing, as well as soil organic carbon.

source
ClimaLand.make_update_driversFunction
make_update_drivers(::AbstractClimaLandDrivers)

Creates and returns a function which updates the driver variables in the default case of no drivers. More generally, this should return a function which updates the driver fields stored in p.drivers.

source
make_update_drivers(driver_tuple)

Creates and returns a function which updates the forcing variables ("drivers"). If no drivers are being used, driver_tuple is empty, and the update function does nothing.

source
make_update_drivers(a::PrescribedGroundConditions{FT}) where {FT}

Creates and returns a function which updates the driver variables in the case of a PrescribedGroundConditions.

source
make_update_drivers(a::PrescribedAtmosphere{FT}) where {FT}

Creates and returns a function which updates the driver variables in the case of a PrescribedAtmosphere.

source
make_update_drivers(a::PrescribedPrecipitation{FT}) where {FT}

Creates and returns a function which updates the driver variables in the case of a PrescribedPrecipitation.

source
make_update_drivers(r::CoupledRadiativeFluxes{FT}) where {FT}

Creates and returns a function which updates the driver variables in the case of a CoupledRadiativeFluxes.

When r.θs is nothing, the cosine zenith angle not changed, and should be updated by the coupler. This differs from the behavior of PrescribedRadiativeFluxes, where the cosine zenith angle is set to NaN if θs is nothing.

Otherwise, the cosine zenith angle is computed using cos.(r.θs(t, r.start_date)).

source
make_update_drivers(r::PrescribedRadiativeFluxes{FT}) where {FT}

Creates and returns a function which updates the driver variables in the case of a PrescribedRadiativeFluxes.

source
make_update_drivers(d::PrescribedSoilOrganicCarbon{FT}) where {FT}

Creates and returns a function which updates the driver variables in the case of a PrescribedSoilOrganicCarbon.

source