Reference

Contents

Index

UrbanTethysChloris.incoming_longwaveMethod
incoming_longwave(Ta::FT, ea::FT, N::FT) where {FT<:AbstractFloat}

Calculate incoming longwave radiation.

Arguments

  • Ta: Air temperature [°C]
  • ea: Vapor pressure [Pa]
  • N: Cloudiness [0-1] or incoming longwave radiation if N > 1 [W/m²]

Returns

  • Latm: Incoming longwave radiation [W/m²]
source
UrbanTethysChloris.set_sun_variablesMethod
set_sun_variables(datam, deltagmt, lon, lat, t_bef, t_aft)

Calculate solar position variables.

Returns tuple of (hS, zetaS, Tsunrise, Tsunset) where:

  • h_S: solar altitude [rad]
  • zeta_S: Sun's azimuth [rad]
  • T_sunrise: sunrise time [h]
  • T_sunset: sunset time [h]
source

Radiation

UrbanTethysChloris.Radiation.RadiationFluxesType
RadiationFluxes{FT<:AbstractFloat}

Structure representing radiation flux components in an urban canyon.

Fields

  • GroundImp: Radiation for impervious ground surface [W/m²]
  • GroundBare: Radiation for bare ground surface [W/m²]
  • GroundVeg: Radiation for vegetated ground surface [W/m²]
  • Tree: Radiation for tree surface [W/m²]
  • WallSun: Radiation for sunlit wall surface [W/m²]
  • WallShade: Radiation for shaded wall surface [W/m²]
  • TotalGround: Total radiation for all ground surfaces [W/m²]
  • TotalCanyon: Total radiation for entire canyon [W/m²]
source
UrbanTethysChloris.Radiation.create_longwave_radiationMethod
create_longwave_radiation(
    A_i::Vector{FT}, fgveg::FT, fgbare::FT, fgimp::FT,
    h_can::FT, w_can::FT, A_t::Union{FT,Nothing},
    Cimp::Bool, Cbare::Bool, Cveg::Bool
) where {FT<:AbstractFloat}

Helper function to create RadiationFluxes objects

source
UrbanTethysChloris.Radiation.direct_shortwave_surfacesMethod
direct_shortwave_surfaces(
    h_can::FT, w_can::FT, d_tree::FT, h_tree::FT, r_tree::FT,
    theta_Z::FT, theta_n::FT, SWR_dir::FT, LAIt::FT, trees::Bool,
    ParVegTree::NamedTuple
) where {FT<:AbstractFloat}

Calculate direct shortwave radiation received by different urban canyon surfaces.

Arguments

  • h_can: normalized building height [-]
  • w_can: normalized street width [-]
  • d_tree: location of trees in the canyon, tree-wall distance [-]
  • h_tree: height of trees, vertical level at the crown center [-]
  • r_tree: size of the tree crown, crown radius [-]
  • theta_Z: solar zenith angle [rad]
  • theta_n: difference between solar azimuth angle and canyon orientation [rad]
  • SWR_dir: direct shortwave radiation W/m^2 of horizontal surfaces [W/m^2]
  • LAIt: leaf area index of the tree [-]
  • trees: boolean indicating if trees are present
  • ParVegTree: named tuple containing vegetation parameters, including Kopt

Returns

  • SWRdir_g: direct shortwave radiation received by the ground [W/m^2]
  • SWRdir_wsun: direct shortwave radiation received by the sunlit wall [W/m^2]
  • SWRdir_wshd: direct shortwave radiation received by the shaded wall [W/m^2]
  • SWRdir_t: direct shortwave radiation received by the tree [W/m^2]
source
UrbanTethysChloris.Radiation.direct_shortwave_treesMethod
direct_shortwave_trees(h_can::FT, d_tree::FT, h_tree::FT, r_tree::FT, theta_Z::FT, theta_n::FT, SWR_dir::FT) where {FT<:AbstractFloat}

Calculate direct shortwave radiation received by two trees in an urban canyon.

Arguments

  • h_can: normalized building height [-]
  • d_tree: location of trees in the canyon, tree-wall distance [-]
  • h_tree: height of trees, vertical level at the crown center [-]
  • r_tree: size of the tree crown, crown radius [-]
  • theta_Z: solar zenith angle [rad]
  • theta_n: difference between solar azimuth angle and canyon orientation [rad]
  • SWR_dir: direct shortwave radiation W/m^2 of horizontal surfaces [W/m^2]

Returns

  • SWR_tree1: direct shortwave radiation received by tree 1 per m^2 tree surface [W/m^2]
  • SWR_tree2: direct shortwave radiation received by tree 2 per m^2 tree surface [W/m^2]
source
UrbanTethysChloris.Radiation.interpolateMethod
interpolate(tree::RadiationFluxes{FT}, notree::RadiationFluxes{FT}, tree_fraction::FT) where {FT<:AbstractFloat}

Combines radiation components from scenarios with and without trees based on the tree fraction.

Arguments

  • tree: RadiationFluxes instance representing scenario with trees
  • notree: RadiationFluxes instance representing scenario without trees
  • tree_fraction: Fraction of area covered by trees [0-1]

Returns

  • RadiationFluxes{FT}: Combined radiation components
source
UrbanTethysChloris.Radiation.longwave_absorbed_no_treeMethod
longwave_absorbed_no_tree(
        h_can::FT, w_can::FT, LWR::FT,
        fgveg::FT, fgbare::FT, fgimp::FT,
        ew::FT, egveg::FT, egbare::FT, egimp::FT,
        Tgimp::FT, Tgbare::FT, Tgveg::FT,
        Twsun::FT, Twshade::FT,
        view_factor::ViewFactor{FT}
) where {FT<:AbstractFloat}

Calculate longwave radiation exchange in an urban canyon without trees.

Arguments

  • h_can: Building height [-]
  • w_can: Ground width [-]
  • LWR: Atmospheric longwave radiation [W/m²]
  • fgveg: Partitioning ground vegetation [-]
  • fgbare: Partitioning ground bare [-]
  • fgimp: Partitioning ground impervious [-]
  • ew: Wall emissivity [-]
  • egveg: Ground vegetation emissivity [-]
  • egbare: Ground bare emissivity [-]
  • egimp: Ground impervious emissivity [-]
  • Tgimp: Temperature of ground impervious [K]
  • Tgbare: Temperature of ground bare [K]
  • Tgveg: Temperature of ground vegetated [K]
  • Twsun: Temperature of wall sun [K]
  • Twshade: Temperature of wall shade [K]
  • view_factor: View factors for radiation exchange

Returns

  • Tuple of (LWRinnT, LWRoutnT, LWRabs_nT) containing longwave radiation components
source
UrbanTethysChloris.Radiation.longwave_absorbed_with_treesMethod
longwave_absorbed_with_trees(
        h_can::FT, w_can::FT, r_tree::FT, LWR::FT,
        fgveg::FT, fgbare::FT, fgimp::FT,
        ew::FT, et::FT, egveg::FT, egbare::FT, egimp::FT,
        Tgimp::FT, Tgbare::FT, Tgveg::FT,
        Twsun::FT, Twshade::FT, Ttree::FT,
        view_factor::ViewFactor{FT}
) where {FT<:AbstractFloat}

Calculate longwave radiation exchange in an urban canyon with trees.

Arguments

  • h_can: Building height [-]
  • w_can: Ground width [-]
  • r_tree: Size of the tree crown, crown radius [-]
  • LWR: Atmospheric longwave radiation [W/m²]
  • fgveg: Partitioning ground vegetation [-]
  • fgbare: Partitioning ground bare [-]
  • fgimp: Partitioning ground impervious [-]
  • ew: Wall emissivity [-]
  • et: Tree emissivity [-]
  • egveg: Ground vegetation emissivity [-]
  • egbare: Ground bare emissivity [-]
  • egimp: Ground impervious emissivity [-]
  • Tgimp: Temperature of ground impervious [K]
  • Tgbare: Temperature of ground bare [K]
  • Tgveg: Temperature of ground vegetated [K]
  • Twsun: Temperature of wall sun [K]
  • Twshade: Temperature of wall shade [K]
  • Ttree: Temperature of tree [K]
  • view_factor: View factors for radiation exchange

Returns

  • Tuple of (LWRinT, LWRoutT, LWRabsT, LWREBT) containing longwave radiation components
source
UrbanTethysChloris.Radiation.shadow_length_no_treeMethod
shadow_length_no_tree(h_can::FT, w_can::FT, theta_Z::FT, theta_n::FT) where {FT<:AbstractFloat}

Calculate shadow lengths in an urban canyon without trees.

Arguments

  • h_can: normalized building height [-]
  • w_can: normalized street width [-]
  • theta_Z: solar zenith angle [rad]
  • theta_n: difference between solar azimuth angle and canyon orientation [rad]

Returns

  • X_Shadow: fraction of ground that is shaded [0-1]
  • X_tree: fraction of ground that is shaded by tree [0-1]
  • n_Shadow: fraction of wall that is shaded [0-1]
  • n_tree: fraction of wall that is shaded by tree [0-1]
source
UrbanTethysChloris.Radiation.shadow_length_with_treesMethod
shadow_length_with_trees(h_can::FT, w_can::FT, d_tree::FT, h_tree::FT, r_tree::FT, theta_Z::FT, theta_n::FT) where {FT<:AbstractFloat}

Calculate shadow lengths in an urban canyon with trees.

Arguments

  • h_can: normalized building height [-]
  • w_can: normalized street width [-]
  • d_tree: location of trees in the canyon, tree-wall distance [-]
  • h_tree: height of trees, vertical level at the crown center [-]
  • r_tree: size of the tree crown, crown radius [-]
  • theta_Z: solar zenith angle [rad]
  • theta_n: difference between solar azimuth angle and canyon orientation [rad]

Returns

  • X_Shadow: fraction of ground that is shaded [0-1]
  • X_tree: fraction of ground that is shaded by tree [0-1]
  • n_Shadow: fraction of wall that is shaded [0-1]
  • n_tree: fraction of wall that is shaded by tree [0-1]
source
UrbanTethysChloris.Radiation.shortwave_absorbed_no_treesMethod
shortwave_absorbed_no_trees(
    h_can::FT, w_can::FT,
    fgveg::FT, fgbare::FT, fgimp::FT,
    awraw::FT, agveg::FT, agbare::FT, agimp::FT,
    SWR_dir::FT, SWR_diff::FT,
    theta_z::FT, theta_n::FT,
    view_factor::ViewFactor{FT},
    ParVegTree::VegetationParams{FT},
    ParWindows::WindowParams{FT},
    bem_enabled::Bool
) where {FT<:AbstractFloat}

Calculate shortwave radiation exchange in an urban canyon without trees.

source
UrbanTethysChloris.Radiation.shortwave_absorbed_with_treesMethod
shortwave_absorbed_with_trees(
        h_can::FT, w_can::FT, h_tree::FT, r_tree::FT, d_tree::FT,
        fgveg::FT, fgbare::FT, fgimp::FT,
        awraw::FT, agveg::FT, agbare::FT, agimp::FT, at::FT,
        LAIt::FT, SWR_dir::FT, SWR_diff::FT,
        theta_Z::FT, theta_n::FT,
        view_factor::ViewFactor{FT},
        ParVegTree::VegetationParams{FT},
        ParWindows::WindowParams{FT},
        bem_enabled::Bool
) where {FT<:AbstractFloat}

Calculate shortwave radiation exchange in an urban canyon with trees.

Returns

Tuple of (SWRinT, SWRoutT, SWRabsT, SWRabsDirT, SWRabsDiffT, SWREBT, albedo_canyon)

source
UrbanTethysChloris.Radiation.total_longwave_absorbedMethod
total_longwave_absorbed(
    temperature_c::AbstractMatrix{FT},
    geometry::UrbanGeometryParameters{FT},
    LWR::FT,
    fractions_ground::LocationSpecificSurfaceFractions{FT},
    prop_optical_ground::VegetatedOpticalProperties{FT},
    prop_optical_wall::SimpleOpticalProperties{FT},
    prop_optical_tree::SimpleOpticalProperties{FT},
    view_factor::ViewFactor{FT}
) where {FT<:AbstractFloat}

Calculate total longwave radiation absorption for urban surfaces.

Arguments

  • temperature_c: Matrix of temperatures for each surface [K]
  • geometry: Urban geometry parameters
  • LWR: Atmospheric longwave radiation [W/m²]
  • fractions_ground: Ground surface fractions
  • prop_optical_ground: Ground optical properties
  • prop_optical_wall: Wall optical properties
  • prop_optical_tree: Tree optical properties
  • view_factor: View factors between surfaces

Returns

Tuple of (LWRint, LWRoutt, LWRabst, LWREBt)

source
UrbanTethysChloris.Radiation.total_shortwave_absorbedMethod
total_shortwave_absorbed(
    geometry::UrbanGeometryParameters{FT},
    SWR_dir::FT,
    SWR_diff::FT,
    theta_n::FT,
    theta_Z::FT,
    fractions_ground::LocationSpecificSurfaceFractions{FT},
    prop_optical_ground::VegetatedOpticalProperties{FT},
    prop_optical_wall::SimpleOpticalProperties{FT},
    prop_optical_tree::SimpleOpticalProperties{FT},
    par_veg_tree::HeightDependentVegetationParameters{FT},
    view_factor::ViewFactor{FT},
    par_windows::WindowParameters{FT},
    bem_enabled::Bool
) where {FT<:AbstractFloat}

Calculate total shortwave radiation absorption for urban surfaces.

Arguments

  • geometry: Urban geometry parameters
  • SWR_dir: Direct shortwave radiation [W/m²]
  • SWR_diff: Diffuse shortwave radiation [W/m²]
  • theta_n: Solar zenith angle [rad]
  • theta_Z: Solar azimuth angle [rad]
  • fractions_ground: Ground surface fractions
  • prop_optical_ground: Ground optical properties
  • prop_optical_wall: Wall optical properties
  • prop_optical_tree: Tree optical properties
  • par_veg_tree: Tree vegetation parameters
  • view_factor: View factors between surfaces
  • par_windows: Window parameters
  • bem_enabled: Whether building energy model is enabled

Returns

Tuple of (SWRint, SWRoutt, SWRabst, SWRabsDirt, SWRabsDifft, SWREBt, albedo_canyon)

source

RayTracing

UrbanTethysChloris.RayTracing.ViewFactorType
ViewFactor{FT<:AbstractFloat} <: AbstractViewFactor{FT}

View factors for radiation exchange between urban surfaces, combining cases with and without trees.

Fields

Without trees (nT)

  • F_gs_nT: Ground to sky view factor
  • F_gw_nT: Ground to wall view factor
  • F_ww_nT: Wall to wall view factor
  • F_wg_nT: Wall to ground view factor
  • F_ws_nT: Wall to sky view factor
  • F_sg_nT: Sky to ground view factor
  • F_sw_nT: Sky to wall view factor

With trees (T)

  • F_gs_T: Ground to sky view factor
  • F_gt_T: Ground to tree view factor
  • F_gw_T: Ground to wall view factor
  • F_ww_T: Wall to wall view factor
  • F_wt_T: Wall to tree view factor
  • F_wg_T: Wall to ground view factor
  • F_ws_T: Wall to sky view factor
  • F_sg_T: Sky to ground view factor
  • F_sw_T: Sky to wall view factor
  • F_st_T: Sky to tree view factor
  • F_tg_T: Tree to ground view factor
  • F_tw_T: Tree to wall view factor
  • F_ts_T: Tree to sky view factor
  • F_tt_T: Tree to tree view factor
source
UrbanTethysChloris.RayTracing.ViewFactorNoTreesType
ViewFactorNoTrees{FT<:AbstractFloat} <: AbstractViewFactor{FT}

View factors for radiation exchange between urban surfaces without trees.

Fields

  • F_gs: Ground to sky view factor
  • F_gw: Ground to wall view factor
  • F_ww: Wall to wall view factor
  • F_wg: Wall to ground view factor
  • F_ws: Wall to sky view factor
  • F_sg: Sky to ground view factor
  • F_sw: Sky to wall view factor
source
UrbanTethysChloris.RayTracing.ViewFactorPointType
ViewFactorPoint{FT<:AbstractFloat} <: AbstractViewFactor{FT}

View factors from a single point to urban surfaces.

Fields

  • F_pg: Point to ground view factor
  • F_ps: Point to sky view factor
  • F_pt: Point to tree view factor
  • F_pwLeft: Point to left wall view factor
  • F_pwRight: Point to right wall view factor
source
UrbanTethysChloris.RayTracing.ViewFactorWithTreesType
ViewFactorWithTrees{FT<:AbstractFloat} <: AbstractViewFactor{FT}

View factors for radiation exchange between urban surfaces with trees.

Fields

  • F_gs: Ground to sky view factor
  • F_gt: Ground to tree view factor
  • F_gw: Ground to wall view factor
  • F_ww: Wall to wall view factor
  • F_wt: Wall to tree view factor
  • F_wg: Wall to ground view factor
  • F_ws: Wall to sky view factor
  • F_sg: Sky to ground view factor
  • F_sw: Sky to wall view factor
  • F_st: Sky to tree view factor
  • F_tg: Tree to ground view factor
  • F_tw: Tree to wall view factor
  • F_ts: Tree to sky view factor
  • F_tt: Tree to tree view factor
source
UrbanTethysChloris.RayTracing.line_segment_intersectMethod
line_segment_intersect(XY1::AbstractMatrix{FT}, XY2::AbstractMatrix{FT})
line_segment_intersect(XY1::SVector{4, FT}, SVector{4, FT})

Find intersections between two sets of line segments in 2D space.

Arguments

  • XY1: vector/matrix where each row is [x1 y1 x2 y2] representing a line segment
  • XY2: vector/matrix where each row is [x1 y1 x2 y2] representing a line segment

Returns

A NamedTuple with fields:

  • intAdjacencyMatrix: N1×N2 boolean matrix indicating intersecting segments
  • intMatrixX: X coordinates of intersection points
  • intMatrixY: Y coordinates of intersection points
  • intNormalizedDistance1To2: Normalized distances from start of XY1 segments to intersections
  • intNormalizedDistance2To1: Normalized distances from start of XY2 segments to intersections
  • parAdjacencyMatrix: Boolean matrix indicating parallel segments
  • coincAdjacencyMatrix: Boolean matrix indicating coincident segments
source
UrbanTethysChloris.RayTracing.view_factors_analyticalMethod
view_factors_analytical(H::FT, W::FT) where {FT<:AbstractFloat}

Compute analytical view factors without trees according to Harman et al. 2004

Arguments

  • H: canyon height [m]
  • W: canyon width [m]

Returns

  • ViewFactorNoTrees: struct containing view factors
source
UrbanTethysChloris.RayTracing.view_factors_canyonMethod
view_factors_canyon(
    geometry::ModelComponents.Parameters.UrbanGeometryParameters{FT},
    person::ModelComponents.Parameters.PersonParameters{FT},
    mc_sample_size::Int = 1000,
    n_rays::Int = 200
) where {FT<:AbstractFloat}

Calculate view factors for an urban street canyon, considering both cases with and without trees.

Arguments

  • geometry: Urban geometry parameters including canyon dimensions and tree properties
  • person: Parameters describing a person's position in the canyon
  • mc_sample_size: Number of Monte Carlo samples for ray tracing (default: 1000)
  • n_rays: Number of rays to emit per sample point (default: 200)

Returns

A tuple containing:

  • vf: ViewFactor object containing both tree and no-tree view factors
  • vf_point: Point-specific view factors
source
UrbanTethysChloris.RayTracing.view_factors_computationMethod
view_factors_computation(
    XSv::Vector{FT}, YSv::Vector{FT}, dmax::FT, sz::FT, dthe::Vector{FT},
    x2::Vector{FT}, z2::Vector{FT}, x3::Vector{FT}, z3::Vector{FT},
    x4::Vector{FT}, z4::Vector{FT}, xc::FT, yc::FT, r::FT, xc2::FT,
    x5::Vector{FT}, z5::Vector{FT}
) where {FT<:AbstractFloat}
view_factors_computation(
    XSv::Vector{FT}, YSv::Vector{FT}, dmax::FT, sz::FT, dthe::Matrix{FT},
    x2::Vector{FT}, z2::Vector{FT}, x3::Vector{FT}, z3::Vector{FT},
    x4::Vector{FT}, z4::Vector{FT}, xc::FT, yc::FT, r::FT, xc2::FT,
    x5::Vector{FT}, z5::Vector{FT}
) where {FT<:AbstractFloat}

Compute view factors for urban canyon surfaces including trees.

Arguments

  • XSv, YSv: Coordinates of emitting points on emitting surface
  • x2,z2: Coordinates of ground surface
  • x3,z3: Coordinates of first wall surface
  • x4,z4: Coordinates of second wall surface
  • xc,yc: Coordinates of centre tree 1
  • xc2: x-coordinate of centre tree 2 (yc2 = yc)
  • r: Radius of tree
  • x5,z5: Coordinates of sky boundary points
  • dmax: Maximum search distance
  • sz: Search step size
  • dthe: Angular resolution in radians

Returns

Tuple of view factors (VG, VW1, VW2, VS, VT1, VT2) for:

  • Ground
  • Wall-1
  • Wall-2
  • Sky
  • Tree-1
  • Tree-2
source
UrbanTethysChloris.RayTracing.view_factors_geometryMethod
view_factors_geometry(H::FT, W::FT, a::FT, ht::FT, d::FT, person::PersonParameters{FT},
                     option_surface::Int, mc_sample_size::Int, n_rays::Int) where {FT<:AbstractFloat}

Compute view factors for urban canyon geometry.

Arguments

  • H: canyon height [m]
  • W: canyon width [m]
  • a: normalized tree radius [-]
  • ht: normalized tree height [-]
  • d: normalized tree distance from wall [-]
  • person: PersonParameters object containing position information
  • option_surface: specifies emitting surface (1=wall1, 2=wall2, 3=ground, 4=tree1, 5=tree2, 6=sky, 7=point)
  • mc_sample_size: number of emitting points per surface
  • n_rays: number of rays emitted per point

Returns

Tuple of view factors (VG, VW1, VW2, VS, VT1, VT2) for ground, walls, sky and trees

source
UrbanTethysChloris.RayTracing.view_factors_ray_tracingMethod
view_factors_ray_tracing(
    H::FT,
    W::FT,
    a::FT,
    ht::FT,
    d::FT,
    person::ModelComponents.Parameters.PersonParameters{FT},
    mc_sample_size::Int,
    n_rays::Int,
) where {FT<:AbstractFloat}

Calculate view factors for urban canyon surfaces using ray tracing method.

Arguments

  • H: canyon height [m]
  • W: canyon width [m]
  • a: normalized tree radius [-]
  • ht: normalized tree height [-]
  • d: normalized tree distance from wall [-]
  • person: PersonParameters object containing position information
  • mc_sample_size: number of Monte Carlo samples for ray tracing
  • n_rays: number of rays to emit per sample point

Returns

Tuple containing:

  • vf: ViewFactorWithTrees object containing view factors between surfaces
  • vfp: ViewFactorPoint object containing point-specific view factors
source
UrbanTethysChloris.RayTracing.view_factors_ray_tracing_reciprocityMethod
view_factors_ray_tracing_reciprocity(H::FT, W::FT, a::FT, ht::FT, d::FT, person::PersonParameters{FT},
                       mc_sample_size::Int, n_rays::Int) where {FT<:AbstractFloat}

Compute reciprocal view factors for urban canyon geometry.

Arguments

  • H: canyon height [m]
  • W: canyon width [m]
  • a: normalized tree radius [-]
  • ht: normalized tree height [-]
  • d: normalized tree distance from wall [-]
  • person: PersonParameters object containing position information
  • mc_sample_size: number of emitting points per surface
  • n_rays: number of rays emitted per point

Returns

Tuple containing (ViewFactorWithTrees, ViewFactorPoint, ViewFactorWithTrees)

source

Soil

UrbanTethysChloris.Soil.conductivity_suctionMethod
conductivity_suction(
    SPAR::Int,
    Ks::FT,
    Osat::FT,
    Ohy::FT,
    L::FT,
    Pe::FT,
    O33::FT,
    alpVG::FT,
    nVG::FT,
    O::FT
) where {FT<:AbstractFloat}

conductivity_suction(
    SPAR::Int,
    Ks::FT,
    Osat::FT,
    Ohy::FT,
    L::FT,
    Pe::FT,
    O33::FT,
    alpVG::FT,
    nVG::FT,
    O::Vector{FT},
) where {FT<:AbstractFloat}

conductivity_suction(
    SPAR::Int,
    Ks::Vector{FT},
    Osat::Vector{FT},
    Ohy::Vector{FT},
    L::Vector{FT},
    Pe::Vector{FT},
    O33::Vector{FT},
    alpVG::Vector{FT},
    nVG::Vector{FT},
    O::Vector{FT},
) where {FT<:AbstractFloat}

weighted_conductivity_suction(
    SPAR::Int,
    Ks::Vector{FT},
    Osat::Vector{FT},
    Ohy::Vector{FT},
    L::Vector{FT},
    Pe::Vector{FT},
    O33::Vector{FT},
    alpVG::Vector{FT},
    nVG::Vector{FT},
    O::FT,
    weights::Vector{FT},
) where {FT<:AbstractFloat}

Calculate soil hydraulic conductivity and matric suction using specified model.

Arguments

  • SPAR: Soil parameterization choice:
    1. Van-Genuchten (1980) corrected
    2. Saxton-Rawls (1986)
  • Ks: Saturated hydraulic conductivity [mm/h]
  • Osat: Saturated water content [m³/m³]
  • Ohy: Hygroscopic water content [m³/m³]
  • L: Lambda parameter for Saxton-Rawls model [-]
  • Pe: Air entry pressure [kPa]
  • O33: Water content at -33 kPa [m³/m³]
  • alpVG: Van Genuchten α parameter [1/mm]
  • nVG: Van Genuchten n parameter [-]
  • O: Current water content [m³/m³]
  • weights: Weights for each soil layer [dimensionless]

Returns

  • Ko: Unsaturated hydraulic conductivity [mm/h]
  • Po: Matric potential [mm]
source
UrbanTethysChloris.Soil.evaporation_layersMethod
evaporation_layers(Zs::Vector{FT}, Zdes::FT)

Calculate the evaporation layer fraction for each soil layer based on the desorption depth.

Arguments

  • Zs::Vector{FT}: Vector of depth layers in mm, from top to bottom [1....m+1]
  • Zdes::FT: Depth of desorption in mm

Returns

  • EvL_Zs::Vector{FT}: Vector of evaporation layer fractions [1...m], where each element represents the fraction of evaporation for the corresponding layer

Notes

  • The function checks if the desorption depth is within the valid range defined by the layers
  • If Zdes < first layer depth: returns zeros with error message "ERROR FIRST LAYER TOO DEPTH"
  • If Zdes > last layer depth: returns zeros with error message "ERROR LAST LAYER TOO SHALLOW"
  • The sum of all fractions equals 1 when Zdes is within valid range
source
UrbanTethysChloris.Soil.infiltrationMethod
infiltration(
    Osat::FT,
    Ohy::FT,
    L::FT,
    alpVG::FT,
    nVG::FT,
    Pe::FT,
    Ks_Zs::FT,
    O33::FT,
    SPAR::Int,
    O::FT,
    Dz::FT,
    WIS::FT,
    cosalp::FT,
    Pond::FT
) where {FT<:AbstractFloat}

Calculate infiltration rates based on soil parameters and conditions.

Arguments

  • Osat: Saturated water content [m³/m³]
  • Ohy: Hygroscopic water content [m³/m³]
  • L: Pore size distribution index [-]
  • alpVG: van Genuchten α parameter [1/mm]
  • nVG: van Genuchten n parameter [-]
  • Pe: Air entry pressure [kPa]
  • Ks_Zs: Saturated hydraulic conductivity [mm/h]
  • O33: Water content at -33 kPa [m³/m³]
  • SPAR: Soil parameterization choice:
    1. van Genuchten (1980)
    2. Saxton-Rawls (1986)
  • O: Current soil water content [m³/m³]
  • Dz: Distance from surface to half-layer [mm]
  • WIS: Water incoming to soil layer [mm/h]
  • cosalp: Cosine of slope angle [-]
  • Pond: Ponding depth [mm]

Returns

  • f: Actual infiltration rate [mm/h]
  • fpot: Potential infiltration rate [mm/h]
source
UrbanTethysChloris.Soil.leakage_bottomMethod
leakage_bottom(
    O::Vector{FT},
    Ks_Zs::Vector{FT},
    Osat::Vector{FT},
    Ohy::Vector{FT},
    L::Vector{FT},
    nVG::Vector{FT},
    Kbot::FT,
    ms::Int,
    SPAR::Int
) where {FT<:AbstractFloat}

Calculate bottom leakage for soil water movement.

Arguments

  • O: Soil water content [m³/m³]
  • Ks_Zs: Saturated hydraulic conductivity [mm/h]
  • Osat: Saturated water content [m³/m³]
  • Ohy: Hygroscopic water content [m³/m³]
  • L: Pore size distribution index [-]
  • nVG: van Genuchten n parameter [-]
  • Kbot: Bottom boundary hydraulic conductivity [mm/h]
  • ms: Index of bottom soil layer [-]
  • SPAR: Soil parameterization choice:
    1. van Genuchten (1980)
    2. Saxton-Rawls (1986)

Returns

  • Lk: Bottom boundary leakage flux [mm/h]
source
UrbanTethysChloris.Soil.root_fraction_generalMethod
root_fraction_general(
    Zs::Vector{FT},
    CASE_ROOT::Int,
    ZR95_H::Vector{FT},
    ZR50_H::Vector{FT},
    ZR95_L::Vector{FT},
    ZR50_L::Vector{FT},
    ZRmax_H::Vector{FT},
    ZRmax_L::Vector{FT}
) where {FT<:AbstractFloat}

Compute root fraction distributions for high and low vegetation using different root profile models.

Arguments

  • Zs: Soil layer depths [mm]
  • CASE_ROOT: Root profile type:
    1. Exponential (Arora and Boer 2005)
    2. Linear Dose Response (Schenk and Jackson 2002)
    3. Constant Profile
    4. Deep (Tap) Root Profile
  • ZR95_H: 95th percentile root depth for high vegetation [mm]
  • ZR50_H: 50th percentile root depth for high vegetation [mm]
  • ZR95_L: 95th percentile root depth for low vegetation [mm]
  • ZR50_L: 50th percentile root depth for low vegetation [mm]
  • ZRmax_H: Maximum root depth for high vegetation [mm]
  • ZRmax_L: Maximum root depth for low vegetation [mm]

Returns

  • RfH_Zs: Root fraction distribution for high vegetation [dimensionless]
  • RfL_Zs: Root fraction distribution for low vegetation [dimensionless]
source
UrbanTethysChloris.Soil.root_soil_conductanceMethod
root_soil_conductance(
    Ks::FT,
    Rl::FT,
    rcyl::FT,
    rroot::FT,
    Zr::FT
) where {FT<:AbstractFloat}

Calculate soil-to-root conductance using the Holtta/Sperry model.

Arguments

  • Ks: Saturated hydraulic conductivity [mm/h]
  • Rl: Root length index [m root/m² ground]
  • rcyl: Radius of the cylinder of soil to which root has access [m]
  • rroot: Root radius [m]
  • Zr: Rooting depth [m]

Returns

  • Ksr: Soil-to-root conductance [mmol H₂O / m² ground s MPa]
source
UrbanTethysChloris.Soil.soil_moisture_conductivity_updateMethod
soil_moisture_conductivity_update(
    V::Vector{FT},
    Pcla::FT,
    Psan::FT,
    Porg::FT,
    Kfc::FT,
    Phy::FT,
    SPAR::Int,
    Kbot::FT,
    CASE_ROOT_H::Int,
    CASE_ROOT_L::Int,
    ZR95_H::Vector{FT},
    ZR95_L::Vector{FT},
    ZR50_H::Vector{FT},
    ZR50_L::Vector{FT},
    ZRmax_H::Vector{FT},
    ZRmax_L::Vector{FT},
    Zs::Vector{FT},
    Rrootl_H::Vector{FT},
    Rrootl_L::Vector{FT},
    PsiL50_H::Vector{FT},
    PsiL50_L::Vector{FT},
    PsiX50_H::Vector{FT},
    PsiX50_L::Vector{FT}
) where {FT<:AbstractFloat}

Updates soil moisture content and calculates hydraulic properties based on current conditions.

Arguments

  • V: Soil water volume per unit area [mm]
  • Pcla: Clay fraction in soil [-]
  • Psan: Sand fraction in soil [-]
  • Porg: Organic matter fraction in soil [-]
  • Kfc: Hydraulic conductivity at field capacity [mm/h]
  • Phy: Soil water potential at hygroscopic point [kPa]
  • SPAR: Soil parameterization choice:
    1. Van-Genuchten (1980)
    2. Saxton-Rawls (1986)
  • Kbot: Hydraulic conductivity at bottom boundary [mm/h]
  • CASE_ROOT_H: Root distribution type for high vegetation [-]
  • CASE_ROOT_L: Root distribution type for low vegetation [-]
  • ZR95_H: 95th percentile root depth for high vegetation [mm]
  • ZR95_L: 95th percentile root depth for low vegetation [mm]
  • ZR50_H: 50th percentile root depth for high vegetation [mm]
  • ZR50_L: 50th percentile root depth for low vegetation [mm]
  • ZRmax_H: Maximum root depth for high vegetation [mm]
  • ZRmax_L: Maximum root depth for low vegetation [mm]
  • Zs: Soil layer depths [mm]
  • Rrootl_H: Root length density for high vegetation [mm/mm³]
  • Rrootl_L: Root length density for low vegetation [mm/mm³]
  • PsiL50_H: Leaf water potential at 50% conductance loss for high vegetation [MPa]
  • PsiL50_L: Leaf water potential at 50% conductance loss for low vegetation [MPa]
  • PsiX50_H: Xylem water potential at 50% conductance loss for high vegetation [MPa]
  • PsiX50_L: Xylem water potential at 50% conductance loss for low vegetation [MPa]

Returns

  • V: Updated soil water volume [mm]
  • O: Volumetric soil moisture content [m³/m³]
  • OS: Surface soil moisture [m³/m³]
  • Psi_soil: Soil water potential [MPa]
  • Psi_s_H: Soil water potential for high vegetation [MPa]
  • Psi_s_L: Soil water potential for low vegetation [MPa]
  • Exwat_H: Water extraction rate for high vegetation [mm/h]
  • Exwat_L: Water extraction rate for low vegetation [mm/h]
  • Ko: Hydraulic conductivity [mm/h]
source
UrbanTethysChloris.Soil.soil_moistures_rich_compMethod
soil_moistures_rich_comp(
    V,
    Lk,
    f,
    EG,
    T_H,
    T_L,
    Qi_in,
    IS,
    SPAR,
    Osat,
    Ohy,
    O33,
    dz,
    Ks_Zs,
    Dz,
    numn,
    L,
    Pe,
    aR,
    aT,
    alpVG,
    nVG,
    cosalp,
    sinalp,
    SN)

Calculate soil moisture changes using Richards equation-based hydrological model.

Arguments

  • V::Vector{FT}: Current soil water volume per unit area [mm]
  • Lk::FT: Bottom leakage rate [mm/h]
  • f::FT: Infiltration rate at top boundary [mm/h]
  • EG::Vector{FT}: Evaporation rates from ground [mm/h]
  • T_H::Vector{FT}: Transpiration rates from hydrophilic roots [mm/h]
  • T_L::Vector{FT}: Transpiration rates from hydrophobic roots [mm/h]
  • Qi_in::Vector{FT}: Lateral inflow rates [mm/h]
  • IS::Vector{FT}: Surface water indicator (1 for surface water present)
  • SPAR::Int: Soil parameter set selector (1: Van Genuchten, 2: Brooks-Corey)
  • Osat::Vector{FT}: Saturated water content at 0 kPa
  • Ohy::Vector{FT}: Residual (hygroscopic) water content
  • O33::Vector{FT}: Water content at -33 kPa tension
  • dz::Vector{FT}: Layer thicknesses [mm]
  • Ks_Zs::Vector{FT}: Saturated hydraulic conductivity [mm/h]
  • Dz::Vector{FT}: Distance between layer centers [mm]
  • numn::Int: Number of soil layers
  • L::Vector{FT}: Lambda parameter (slope of log tension-moisture curve)
  • Pe::Vector{FT}: Air entry tension (bubbling pressure) [kPa]
  • aR::FT: Horizontal length scale parameter
  • aT::FT: Terrain curvature parameter
  • alpVG::Vector{FT}: van Genuchten α parameter
  • nVG::Vector{FT}: van Genuchten n parameter
  • cosalp::FT: Cosine of slope angle
  • sinalp::FT: Sine of slope angle
  • SN::FT: Stream network identifier

Returns

  • dV::Vector{FT}: Change rates of soil water volume per unit area [mm/h]
source
UrbanTethysChloris.Soil.soil_moistures_rich_comp_lat2Method
soil_moistures_rich_comp_lat2(
    Vlat::Vector{FT},
    dz::FT,
    SPAR::Int,
    Ks::FT,
    Osat::FT,
    Ohy::FT,
    L::FT,
    Pe::FT,
    O33::FT,
    alpVG::FT,
    nVG::FT,
    C1::FT,
    C2::FT,
    f1::FT,
    f2::FT,
    Wcan::FT
) where {FT<:AbstractFloat}

Calculate lateral soil moisture redistribution between soil columns.

Arguments

  • Vlat: Vector of soil water volumes per unit area for both columns [mm]
  • dz: Layer thickness [mm]
  • SPAR: Soil parameterization choice (1: Van-Genuchten, 2: Saxton-Rawls)
  • Ks: Saturated hydraulic conductivity [mm/h]
  • Osat: Saturated water content [m³/m³]
  • Ohy: Hygroscopic water content [m³/m³]
  • L: Lambda parameter (for Saxton-Rawls model) [-]
  • Pe: Air entry pressure [kPa]
  • O33: Water content at -33 kPa [m³/m³]
  • alpVG: Van Genuchten α parameter [1/mm]
  • nVG: Van Genuchten n parameter [-]
  • C1: Column 1 coefficient [-]
  • C2: Column 2 coefficient [-]
  • f1: Fraction for column 1 [-]
  • f2: Fraction for column 2 [-]
  • Wcan: Canopy width [m]

Returns

  • dVlat: Change rates of soil water volume per unit area for both columns [mm/h]

Notes

  • If f1 or f2 is zero, the corresponding flux is set to zero.
source
UrbanTethysChloris.Soil.soil_moistures_rich_comp_lat3Method
soil_moistures_rich_comp_lat3(
    Vlat::Vector{FT},
    dz::FT,
    SPAR::Int,
    Ks::FT,
    Osat::FT,
    Ohy::FT,
    L::FT,
    Pe::FT,
    O33::FT,
    alpVG::FT,
    nVG::FT,
    Cimp::FT,
    Cbare::FT,
    Cveg::FT,
    fimp::FT,
    fbare::FT,
    fveg::FT,
    Wcan::FT
) where {FT<:AbstractFloat}

Calculate lateral soil moisture redistribution between three soil columns.

Arguments

  • Vlat: Soil water volume per unit area for each column [mm]
  • dz: Soil layer thickness [mm]
  • SPAR: Soil parameterization choice:
    1. van Genuchten (1980)
    2. Saxton-Rawls (1986)
  • Ks: Saturated hydraulic conductivity [mm/h]
  • Osat: Saturated water content [m³/m³]
  • Ohy: Hygroscopic water content [m³/m³]
  • L: Pore size distribution index [-]
  • Pe: Air entry pressure [kPa]
  • O33: Water content at -33 kPa [m³/m³]
  • alpVG: van Genuchten α parameter [1/mm]
  • nVG: van Genuchten n parameter [-]
  • Cimp: Contact length of impervious column [mm]
  • Cbare: Contact length of bare soil column [mm]
  • Cveg: Contact length of vegetated column [mm]
  • fimp: Fraction of impervious area [-]
  • fbare: Fraction of bare soil area [-]
  • fveg: Fraction of vegetated area [-]
  • Wcan: Canyon width [mm]

Returns

  • dVlat::Vector{FT}: Change rates of soil water volume [mm/h] for:
    1. Impervious column
    2. Bare soil column
    3. Vegetated column
source
UrbanTethysChloris.Soil.soil_parametersMethod
soil_parameters(
    Psan::FT,
    Pcla::FT,
    Porg::FT
) where {FT<:AbstractFloat}

Calculate soil hydraulic and thermal parameters based on soil texture composition.

Arguments

  • Psan: Sand fraction in soil [0-1]
  • Pcla: Clay fraction in soil [0-1]
  • Porg: Organic matter fraction in soil [0-1]

Returns

  • Osat::FT: Saturated soil water content [-]
  • L::FT: Slope of logarithmic tension-moisture curve [-]
  • Pe::FT: Air entry tension [kPa]
  • Ks::FT: Saturated hydraulic conductivity [mm/h]
  • O33::FT: Soil water content at 33 kPa [-]
  • rsd::FT: Soil bulk density [kg/m³]
  • lan_dry::FT: Thermal conductivity of dry soil [W/m K]
  • lan_s::FT: Thermal conductivity of solid soil components [W/m K]
  • cv_s::FT: Volumetric heat capacity of solid soil components [J/m³ K]
  • K_usle::FT: USLE K factor for soil erodibility [kgh/Jmm]
source
UrbanTethysChloris.Soil.soil_parameters2Method
soil_parameters2(
    Osat::Vector{FT},
    L::Vector{FT},
    Pe::Vector{FT},
    Ks::Vector{FT},
    O33::Vector{FT},
    nVG::Vector{FT},
    alpVG::Vector{FT},
    Kfc::FT,
    Pss::FT,
    Pwp::FT,
    Phy::FT,
    Ohy::Vector{FT}=zeros(FT, length(Osat)),
    SPAR::Int=2
) where {FT<:AbstractFloat}

Calculate soil moisture retention points based on soil hydraulic parameters.

Arguments

  • Osat::Vector{FT}: Saturated soil water content [-]
  • L::Vector{FT}: Pore size distribution index [-]
  • Pe::Vector{FT}: Air entry pressure [kPa]
  • Ks::Vector{FT}: Saturated hydraulic conductivity [mm/h]
  • O33::Vector{FT}: Soil water content at -33 kPa [-]
  • nVG::Vector{FT}: van Genuchten n parameter [-]
  • alpVG::Vector{FT}: van Genuchten α parameter [1/mm]
  • Kfc::FT: Hydraulic conductivity at field capacity [mm/h]
  • Pss::FT: Soil water potential at steady state [kPa]
  • Pwp::FT: Permanent wilting point potential [kPa]
  • Phy::FT: Hygroscopic point potential [kPa]
  • Ohy::Vector{FT}: Initial hygroscopic water content [-]
  • SPAR::Int: Soil parameter set:
    1. van Genuchten (1980)
    2. Brooks-Corey (default)

Returns

  • Ofc::Vector{FT}: Soil water content at field capacity [-]
  • Oss::Vector{FT}: Soil water content at steady state [-]
  • Owp::Vector{FT}: Soil water content at wilting point [-]
  • Ohy::Vector{FT}: Soil water content at hygroscopic point [-]
source
UrbanTethysChloris.Soil.soil_parameters_totalMethod
soil_parameters_total(
    Pcla::FT,
    Psan::FT,
    Porg::FT,
    Kfc::FT,
    Phy::FT,
    SPAR::Int,
    Kbot::FT,
    CASE_ROOT_H::Int,
    CASE_ROOT_L::Int,
    ZR95_H::Vector{FT},
    ZR95_L::Vector{FT},
    ZR50_H::Vector{FT},
    ZR50_L::Vector{FT},
    ZRmax_H::Vector{FT},
    ZRmax_L::Vector{FT},
    Zs::Vector{FT}
) where {FT<:AbstractFloat}

Calculate total soil parameters based on composition and root distribution.

Arguments

  • Pcla: Clay fraction in soil [-]
  • Psan: Sand fraction in soil [-]
  • Porg: Organic matter fraction in soil [-]
  • Kfc: Hydraulic conductivity at field capacity [mm/h]
  • Phy: Soil water potential at hygroscopic point [kPa]
  • SPAR: Soil parameterization choice:
    1. van Genuchten (1980)
    2. Saxton-Rawls (1986)
  • Kbot: Hydraulic conductivity at bottom boundary [mm/h]
  • CASE_ROOT_H: Root distribution type for high vegetation [-]
  • CASE_ROOT_L: Root distribution type for low vegetation [-]
  • ZR95_H: 95th percentile root depth for high vegetation [mm]
  • ZR95_L: 95th percentile root depth for low vegetation [mm]
  • ZR50_H: 50th percentile root depth for high vegetation [mm]
  • ZR50_L: 50th percentile root depth for low vegetation [mm]
  • ZRmax_H: Maximum root depth for high vegetation [mm]
  • ZRmax_L: Maximum root depth for low vegetation [mm]
  • Zs: Soil layer depths [mm]

Returns

  • Zs: Soil layer depths [mm]
  • dz: Layer thicknesses [mm]
  • ms: Number of soil layers [-]
  • Osat: Saturated water content [m³/m³]
  • Ohy: Hygroscopic water content [m³/m³]
  • nVG: van Genuchten n parameter [-]
  • alpVG: van Genuchten α parameter [1/mm]
  • Ks_Zs: Saturated hydraulic conductivity [mm/h]
  • L: Pore size distribution index [-]
  • Pe: Air entry pressure [kPa]
  • O33: Water content at -33 kPa [m³/m³]
  • SPAR: Soil parameterization type [-]
  • EvL_Zs: Evaporation layer fractions [-]
  • Inf_Zs: Infiltration layer fractions [-]
  • RfH_Zs: Root fraction distribution for high vegetation [-]
  • RfL_Zs: Root fraction distribution for low vegetation [-]
  • Zinf: Infiltration depth [mm]
  • Kbot: Bottom boundary conductivity [mm/h]
  • Slo_pot: Slope fractions [-]
  • Dz: Layer center distances [mm]
  • aR: Horizontal length scale [-]
  • aTop: Area to contour length ratio [mm]
  • rsd: Dry soil density [kg/m³]
  • lan_dry: Dry soil thermal conductivity [W/m K]
  • lan_s: Solid soil thermal conductivity [W/m K]
  • cv_s: Solid soil volumetric heat capacity [J/m³ K]
source
UrbanTethysChloris.Soil.soil_thermal_propertiesMethod
soil_thermal_properties(Tdp, rsd, lan_dry, lan_s, cv_s, Osat, Ohy, O) where {FT<:AbstractFloat}

Compute soil thermal properties based on soil moisture and composition.

Arguments

  • Tdp::Vector{FT},: Dampening temperature [°C]
  • rsd::Vector{FT}: Density dry soil [kg/m³]
  • lan_dry::Vector{FT}: Thermal conductivity dry soil [W/m K]
  • lan_s::Vector{FT}: Thermal conductivity soil solid [W/m K]
  • cv_s::Vector{FT}: Volumetric heat capacity soil solid [J/m³ K]
  • Osat::Vector{FT}: Water content at saturation [-]
  • Ohy::Vector{FT}: Residual/hygroscopic water content [-]
  • O::Vector{FT}: Soil moisture at previous time step [-]

Returns

  • lanS::Vector{FT}: Thermal conductivity of soil [W/m K]
  • cv_Soil::Vector{FT}: Volumetric heat capacity of soil [J/m³ K]
  • CTt::Vector{FT}: Total thermal capacity of soil [K m²/J]
source
UrbanTethysChloris.Soil.soil_water_multilayerMethod
soil_water_multilayer(
    V::Vector{FT},
    Zs::Vector{FT},
    dz::Vector{FT},
    n::Int,
    Osat::Vector{FT},
    Ohy::Vector{FT},
    nVG::Vector{FT},
    alpVG::Vector{FT},
    Ks_Zs::Vector{FT},
    L::Vector{FT},
    Pe::Vector{FT},
    O33::Vector{FT},
    SPAR::Int,
    EvL_Zs::Vector{FT},
    Inf_Zs::Vector{FT},
    RfH_Zs::Matrix{FT},
    RfL_Zs::Matrix{FT},
    Rrootl_H::Vector{FT},
    Rrootl_L::Vector{FT},
    PsiL50_H::Vector{FT},
    PsiL50_L::Vector{FT},
    PsiX50_H::Vector{FT},
    PsiX50_L::Vector{FT}
) where {FT<:AbstractFloat}

Calculate soil water dynamics and potential water extraction rates.

Arguments

  • V: Soil water volume per unit area [mm]
  • Zs: Soil layer depths [mm]
  • dz: Layer thicknesses [mm]
  • n: Number of soil layers [-]
  • Osat: Saturated water content [m³/m³]
  • Ohy: Hygroscopic water content [m³/m³]
  • nVG: van Genuchten n parameter [-]
  • alpVG: van Genuchten α parameter [1/mm]
  • Ks_Zs: Saturated hydraulic conductivity [mm/h]
  • L: Pore size distribution index [-]
  • Pe: Air entry pressure [kPa]
  • O33: Water content at -33 kPa [m³/m³]
  • SPAR: Soil parameterization choice:
    1. van Genuchten (1980)
    2. Saxton-Rawls (1986)
  • EvL_Zs: Evaporation layer fractions [-]
  • Inf_Zs: Infiltration layer fractions [-]
  • RfH_Zs: Root fraction distribution for high vegetation [-]
  • RfL_Zs: Root fraction distribution for low vegetation [-]
  • Rrootl_H: Root length density for high vegetation [mm/mm³]
  • Rrootl_L: Root length density for low vegetation [mm/mm³]
  • PsiL50_H: Leaf water potential at 50% loss for high vegetation [MPa]
  • PsiL50_L: Leaf water potential at 50% loss for low vegetation [MPa]
  • PsiX50_H: Xylem water potential at 50% loss for high vegetation [MPa]
  • PsiX50_L: Xylem water potential at 50% loss for low vegetation [MPa]

Returns

  • O: Volumetric soil moisture content [m³/m³]
  • ZWT: Water table depth [mm]
  • OF: Infiltration water content [-]
  • OS: Surface soil moisture [-]
  • Psi_s_H: Soil water potential for high vegetation [MPa]
  • Psi_s_L: Soil water potential for low vegetation [MPa]
  • gsr_H: Root-soil conductance for high vegetation [mmol/m²/s/MPa]
  • gsr_L: Root-soil conductance for low vegetation [mmol/m²/s/MPa]
  • Exwat_H: Water extraction rate for high vegetation [mm/h]
  • Exwat_L: Water extraction rate for low vegetation [mm/h]
  • Rd: Surface runoff [mm]
  • WTR: Water table rise [mm]
  • POT: Soil water potential [mm]
  • OH: Average soil moisture for high vegetation [-]
  • OL: Average soil moisture for low vegetation [-]
source
UrbanTethysChloris.Soil.volume_correctionMethod
volume_correction(V, EvL_Zs, RfH_Zs, RfL_Zs, EG, T_H, T_L, Lk)

Volume correction for negative values.

Arguments

  • V: Volume array
  • EvL_Zs: Evaporation layers array
  • RfH_Zs: Root fraction high vegetation layers array
  • RfL_Zs: Root fraction low vegetation layers array
  • EG: Ground evaporation
  • T_H: High vegetation transpiration array
  • T_L: Low vegetation transpiration array
  • Lk: Leakage value

Returns

  • V: Corrected volume array
  • T_H: Updated high vegetation transpiration
  • T_L: Updated low vegetation transpiration
  • EG: Updated ground evaporation
  • Lk: Updated leakage
source

Water

UrbanTethysChloris.Water.water_canyonMethod
water_canyon(
    MeteoData::NamedTuple,
    Int_ittm::NamedTuple,
    Owater_ittm::NamedTuple,
    Runon_ittm::NamedTuple,
    Qinlat_ittm::NamedTuple,
    Etree_In::FT,
    Egveg_In::FT,
    Egimp_Pond::FT,
    Egbare_Pond::FT,
    Egveg_Pond::FT,
    Egbare_soil::FT,
    Egveg_soil::FT,
    TEgveg::FT,
    TEtree::FT,
    ParSoilGround::NamedTuple,
    ParInterceptionTree::NamedTuple,
    ParCalculation::NamedTuple,
    ParVegGround::NamedTuple,
    ParVegTree::NamedTuple,
    FractionsGround::NamedTuple,
    geometry::NamedTuple,
    ParTree::NamedTuple,
    Gemeotry_m::NamedTuple,
    Anthropogenic::NamedTuple
) where {FT<:AbstractFloat}

Calculate water balance for canyon surfaces including trees, vegetation, and impervious areas.

Arguments

  • MeteoData: Meteorological data
  • Int_ittm: Previous timestep interception values
  • Owater_ittm: Previous timestep soil moisture values
  • Runon_ittm: Previous timestep runon values
  • Qinlat_ittm: Previous timestep lateral flow values
  • Etree_In: Tree evaporation [kg/m²s]
  • Egveg_In: Ground vegetation evaporation [kg/m²s]
  • Egimp_Pond: Impervious ground evaporation [kg/m²s]
  • Egbare_Pond: Bare ground evaporation [kg/m²s]
  • Egveg_Pond: Vegetated ground evaporation [kg/m²s]
  • Egbare_soil: Bare soil evaporation [kg/m²s]
  • Egveg_soil: Vegetated soil evaporation [kg/m²s]
  • TEgveg: Ground vegetation transpiration [kg/m²s]
  • TEtree: Tree transpiration [kg/m²s]
  • ParSoilGround: Soil parameters for ground
  • ParInterceptionTree: Tree interception parameters
  • ParCalculation: Calculation parameters
  • ParVegGround: Ground vegetation parameters
  • ParVegTree: Tree vegetation parameters
  • FractionsGround: Ground surface fractions
  • geometry: Geometry parameters
  • ParTree: Tree parameters
  • Gemeotry_m: Geometry parameters in meters
  • Anthropogenic: Anthropogenic water inputs

Returns

  • q_tree_dwn: Tree throughfall [mm/dth]
  • In_tree: Tree interception [mm]
  • dIn_tree_dt: Change in tree interception [mm/dth]
  • q_gveg_dwn: Ground vegetation throughfall [mm/dth]
  • In_gveg: Ground vegetation interception [mm]
  • dIn_gveg_dt: Change in ground vegetation interception [mm/dth]
  • q_gimp_runoff: Impervious ground runoff [mm/dth]
  • In_gimp: Impervious ground interception [mm]
  • dIn_gimp_dt: Change in impervious ground interception [mm/dth]
  • f_inf_gimp: Impervious ground infiltration [mm/h]
  • q_gbare_runoff: Bare ground runoff [mm/dth]
  • In_gbare: Bare ground interception [mm]
  • dIn_gbare_dt: Change in bare ground interception [mm/dth]
  • f_inf_gbare: Bare ground infiltration [mm/h]
  • q_gveg_runoff: Vegetated ground runoff [mm/dth]
  • In_gveg_pond: Ground vegetation ponding [mm]
  • dIn_gveg_pond_dt: Change in ground vegetation ponding [mm/dth]
  • f_inf_gveg: Vegetated ground infiltration [mm/h]
  • V_gimp: Impervious ground water volume [mm]
  • O_gimp: Impervious ground water content [-]
  • OS_gimp: Impervious ground surface water content [-]
  • Lk_gimp: Impervious ground leakage [mm/h]
  • Psi_s_H_gimp: Impervious ground high soil water potential [MPa]
  • Psi_s_L_gimp: Impervious ground low soil water potential [MPa]
  • Exwat_H_gimp: Impervious ground high extractable water [mm]
  • Exwat_L_gimp: Impervious ground low extractable water [mm]
  • Rd_gimp: Impervious ground surface runoff [mm]
  • TEgveg_imp: Ground vegetation transpiration from impervious [kg/m²s]
  • TEtree_imp: Tree transpiration from impervious [kg/m²s]
  • Egimp_soil: Impervious soil evaporation [kg/m²s]
  • dV_dt_gimp: Change in impervious ground water volume [mm/dth]
  • Psi_soil_gimp: Impervious ground soil water potential [MPa]
  • Kf_gimp: Impervious ground hydraulic conductivity [mm/h]
  • V_gbare: Bare ground water volume [mm]
  • O_gbare: Bare ground water content [-]
  • OS_gbare: Bare ground surface water content [-]
  • Lk_gbare: Bare ground leakage [mm]
  • Psi_s_H_gbare: Bare ground high soil water potential [MPa]
  • Psi_s_L_gbare: Bare ground low soil water potential [MPa]
  • Exwat_H_gbare: Bare ground high extractable water [mm]
  • Exwat_L_gbare: Bare ground low extractable water [mm]
  • Rd_gbare: Bare ground surface runoff [mm]
  • TEgveg_bare: Ground vegetation transpiration from bare [kg/m²s]
  • TEtree_bare: Tree transpiration from bare [kg/m²s]
  • Egbare_Soil: Bare soil evaporation [kg/m²s]
  • dV_dt_gbare: Change in bare ground water volume [mm/dth]
  • Psi_soil_gbare: Bare ground soil water potential [MPa]
  • Kf_gbare: Bare ground hydraulic conductivity [mm/h]
  • V_gveg: Vegetated ground water volume [mm]
  • O_gveg: Vegetated ground water content [-]
  • OS_gveg: Vegetated ground surface water content [-]
  • Lk_gveg: Vegetated ground leakage [mm]
  • Psi_s_H_gveg: Vegetated ground high soil water potential [MPa]
  • Psi_s_L_gveg: Vegetated ground low soil water potential [MPa]
  • Exwat_H_gveg: Vegetated ground high extractable water [mm]
  • Exwat_L_gveg: Vegetated ground low extractable water [mm]
  • Rd_gveg: Vegetated ground surface runoff [mm]
  • TEgveg_veg: Ground vegetation transpiration from vegetated [kg/m²s]
  • TEtree_veg: Tree transpiration from vegetated [kg/m²s]
  • Egveg_Soil: Vegetated soil evaporation [kg/m²s]
  • dV_dt_gveg: Change in vegetated ground water volume [mm/dth]
  • Psi_soil_gveg: Vegetated ground soil water potential [MPa]
  • Kf_gveg: Vegetated ground hydraulic conductivity [mm/h]
  • Qin_imp: Lateral flow to impervious ground [mm/dth]
  • Qin_bare: Lateral flow to bare ground [mm/dth]
  • Qin_veg: Lateral flow to vegetated ground [mm/dth]
  • Qin_bare2imp: Lateral flow from bare to impervious [mm/dth]
  • Qin_bare2veg: Lateral flow from bare to vegetated [mm/dth]
  • Qin_imp2bare: Lateral flow from impervious to bare [mm/dth]
  • Qin_imp2veg: Lateral flow from impervious to vegetated [mm/dth]
  • Qin_veg2imp: Lateral flow from vegetated to impervious [mm/dth]
  • Qin_veg2bare: Lateral flow from vegetated to bare [mm/dth]
  • V: Total water volume [mm]
  • O: Total water content [-]
  • OS: Surface water content [-]
  • Lk: Total leakage [mm/h]
  • Rd: Surface runoff [mm]
  • dV_dt: Change in total water volume [mm/dth]
  • Psi_s_L: Low soil water potential [MPa]
  • Exwat_L: Low extractable water [mm]
  • TEgveg_tot: Total ground vegetation transpiration [kg/m²s]
  • Psi_s_H_tot: Total high soil water potential [MPa]
  • Exwat_H: Total extractable water [mm]
  • TEtree_tot: Total tree transpiration [kg/m²s]
  • EB_TEtree: Tree transpiration energy balance [W/m²]
  • EB_TEgveg: Ground vegetation transpiration energy balance [W/m²]
  • WBIndv: Individual water balance components
  • WBTot: Total water balance components
  • Runoff: Total runoff [mm/dth]
  • Runon: Total runon [mm/dth]
  • Etot: Total evapotranspiration [mm/dth]
  • DeepGLk: Deep ground leakage [mm/dth]
  • StorageTot: Total water storage change [mm/dth]
source
UrbanTethysChloris.Water.water_groundMethod
water_ground(
    q_runon_veg::FT,
    Runon_tm1::FT,
    E_ground::FT,
    Otm1::Vector{FT},
    In_ground_tm1::FT,
    In_max_ground::FT,
    Pcla::FT,
    Psan::FT,
    Porg::FT,
    Kfc::FT,
    Phy::FT,
    SPAR::Int,
    Kbot::FT,
    CASE_ROOT_H::Int,
    CASE_ROOT_L::Int,
    ZR95_H::Vector{FT},
    ZR95_L::Vector{FT},
    ZR50_H::Vector{FT},
    ZR50_L::Vector{FT},
    ZRmax_H::Vector{FT},
    ZRmax_L::Vector{FT},
    Zs::Vector{FT},
    dth::FT,
    row::FT
) where {FT<:AbstractFloat}

Calculate ground water dynamics and potential infiltration rates.

Arguments

  • q_runon_veg: Runon from vegetation [mm/dth]
  • Runon_tm1: Previous timestep runon [mm/dth]
  • E_ground: Ground evaporation [kg/m²s]
  • Otm1: Previous timestep soil moisture [m³/m³]
  • In_ground_tm1: Previous timestep ground interception [mm]
  • In_max_ground: Maximum ground interception capacity [mm]
  • Pcla: Clay fraction in soil [-]
  • Psan: Sand fraction in soil [-]
  • Porg: Organic matter fraction in soil [-]
  • Kfc: Hydraulic conductivity at field capacity [mm/h]
  • Phy: Soil water potential at hygroscopic point [kPa]
  • SPAR: Soil parameterization choice [-]
  • Kbot: Hydraulic conductivity at bottom boundary [mm/h]
  • CASE_ROOT_H: Root distribution type for high vegetation [-]
  • CASE_ROOT_L: Root distribution type for low vegetation [-]
  • ZR95_H: 95th percentile root depth for high vegetation [mm]
  • ZR95_L: 95th percentile root depth for low vegetation [mm]
  • ZR50_H: 50th percentile root depth for high vegetation [mm]
  • ZR50_L: 50th percentile root depth for low vegetation [mm]
  • ZRmax_H: Maximum root depth for high vegetation [mm]
  • ZRmax_L: Maximum root depth for low vegetation [mm]
  • Zs: Soil layer depths [mm]
  • dth: Calculation time step [h]
  • row: Water density [kg/m³]

Returns

  • q_runon_ground: Ground runoff [mm/dth]
  • In_ground: Ground interception [mm]
  • dIn_ground_dt: Change in ground interception [mm/dth]
  • f_ground: Infiltration rate [mm/h]
  • WBalance_In_ground: Water balance check [mm/dth]
source
UrbanTethysChloris.Water.water_imperviousMethod
water_impervious(
    Rain::FT,
    Runon_tm1::FT,
    E_imp::FT,
    In_imp_tm1::FT,
    dth::FT,
    row::FT,
    In_max_imp::FT,
    K_imp::FT
) where {FT<:AbstractFloat}

Calculates water balance for impervious surfaces.

Arguments

  • Rain: precipitation [mm/dth]
  • Runon_tm1: Previous timestep runon [mm/dth]
  • E_imp: evaporation from impervious surfaces [kg/m²s]
  • In_imp_tm1: Interception from previous time step [mm]
  • dth: calculation time step [h]
  • row: density of water [kg/m³]
  • In_max_imp: Maximum interception capacity of urban area [mm]
  • K_imp: Hydraulic conductivity of impervious area [mm/h]

Returns

  • q_runon_imp: Runoff [mm/dth]
  • In_imp: Interception [mm]
  • dIn_imp_dt: Change in interception [mm/dth]
  • Lk_imp: Leakage from impervious area [mm/h]
  • WBalance_In_imp: Water balance check [mm/dth]
source
UrbanTethysChloris.Water.water_roofMethod
water_roof(
    Eroof_imp::FT,
    Eroof_veg::FT,
    Eroof_ground::FT,
    Eroof_soil::FT,
    TEroof_veg::FT,
    MeteoData::NamedTuple,
    Int_ittm::NamedTuple,
    Owater_ittm::NamedTuple,
    Runon_ittm::NamedTuple,
    FractionsRoof::NamedTuple,
    ParSoilRoof::NamedTuple,
    ParCalculation::NamedTuple,
    ParVegRoof::NamedTuple,
    Anthropogenic::NamedTuple
) where {FT<:AbstractFloat}

Calculate water balance for roof surfaces including vegetation, impervious areas, and soil.

Arguments

  • Eroof_imp: Evaporation from impervious surfaces [kg/m²s]
  • Eroof_veg: Evaporation from vegetation [kg/m²s]
  • Eroof_ground: Evaporation from ground under vegetation [kg/m²s]
  • Eroof_soil: Evaporation from soil [kg/m²s]
  • TEroof_veg: Transpiration from vegetation [kg/m²s]
  • MeteoData: Meteorological data
  • Int_ittm: Previous timestep interception values
  • Owater_ittm: Previous timestep soil moisture values
  • Runon_ittm: Previous timestep runon values
  • FractionsRoof: Roof surface fractions
  • ParSoilRoof: Soil parameters for roof
  • ParCalculation: Calculation parameters
  • ParVegRoof: Vegetation parameters for roof
  • Anthropogenic: Anthropogenic water inputs

Returns

  • q_runon_imp: Runoff from impervious surfaces [mm/dth]
  • In_imp: Interception on impervious surfaces [mm]
  • dIn_imp_dt: Change in impervious surface interception [mm/dth]
  • Lk_imp: Leakage from impervious surfaces [mm/h]
  • q_runon_veg: Runoff from vegetation [mm/dth]
  • In_veg: Vegetation interception [mm]
  • dIn_veg_dt: Change in vegetation interception [mm/dth]
  • q_runon_ground: Ground runoff [mm/dth]
  • In_ground: Ground interception [mm]
  • dIn_ground_dt: Change in ground interception [mm/dth]
  • dV_dt: Change in soil water volume [mm/dth]
  • f_ground: Ground infiltration rate [mm/h]
  • V: Soil water volume [mm]
  • O: Soil moisture [-]
  • OS: Surface soil moisture [-]
  • Lk: Soil leakage [mm/h]
  • Psi_s: Soil water potential [MPa]
  • Exwat: Extractable water [mm/h]
  • Rd: Surface runoff [mm/dth]
  • TEroof_veg: Updated vegetation transpiration [kg/m²s]
  • Eroof_soil: Updated soil evaporation [kg/m²s]
  • Runoff: Total roof runoff [mm/dth]
  • Runon: Total roof runon [mm/dth]
  • WBalance_In_imp: Water balance for impervious surfaces [mm/dth]
  • WBalance_In_veg: Water balance for vegetation [mm/dth]
  • WBalance_In_ground: Water balance for ground [mm/dth]
  • WBalance_soil: Water balance for soil [mm/dth]
  • WBalance_imp_tot: Total water balance for impervious surfaces [mm/dth]
  • WBalance_veg_tot: Total water balance for vegetation [mm/dth]
  • WBalance_tot: Total water balance [mm/dth]
source
UrbanTethysChloris.Water.water_soilMethod
water_soil(
    Otm1::Vector{FT},
    f::FT,
    TE_H::Vector{FT},
    TE_L::Vector{FT},
    E_soil::FT,
    Qlat_in::Vector{FT},
    dth::FT,
    Pcla::FT,
    Psan::FT,
    Porg::FT,
    Kfc::FT,
    Phy::FT,
    SPAR::Int,
    Kbot::FT,
    CASE_ROOT_H::Int,
    CASE_ROOT_L::Int,
    ZR95_H::Vector{FT},
    ZR95_L::Vector{FT},
    ZR50_H::Vector{FT},
    ZR50_L::Vector{FT},
    ZRmax_H::Vector{FT},
    ZRmax_L::Vector{FT},
    Rrootl_H::Vector{FT},
    Rrootl_L::Vector{FT},
    PsiL50_H::Vector{FT},
    PsiL50_L::Vector{FT},
    PsiX50_H::Vector{FT},
    PsiX50_L::Vector{FT},
    Zs::Vector{FT},
    row::FT
) where {FT<:AbstractFloat}

Calculate soil water dynamics and hydrological fluxes.

Arguments

  • Otm1: Previous timestep soil moisture [-]
  • f: Infiltration rate [mm/h]
  • TE_H: High vegetation transpiration [kg/m²s]
  • TE_L: Low vegetation transpiration [kg/m²s]
  • E_soil: Soil evaporation [kg/m²s]
  • Qlat_in: Lateral water flux [mm/dth]
  • dth: Time step [h]
  • Pcla: Clay fraction [-]
  • Psan: Sand fraction [-]
  • Porg: Organic matter fraction [-]
  • Kfc: Hydraulic conductivity at field capacity [mm/h]
  • Phy: Soil water potential at hygroscopic point [kPa]
  • SPAR: Soil parameterization choice [-]
  • Kbot: Bottom boundary conductivity [mm/h]
  • CASE_ROOT_H: High vegetation root distribution type [-]
  • CASE_ROOT_L: Low vegetation root distribution type [-]
  • ZR95_H: 95th percentile root depth for high vegetation [mm]
  • ZR95_L: 95th percentile root depth for low vegetation [mm]
  • ZR50_H: 50th percentile root depth for high vegetation [mm]
  • ZR50_L: 50th percentile root depth for low vegetation [mm]
  • ZRmax_H: Maximum root depth for high vegetation [mm]
  • ZRmax_L: Maximum root depth for low vegetation [mm]
  • Rrootl_H: Root length density for high vegetation [mm/mm³]
  • Rrootl_L: Root length density for low vegetation [mm/mm³]
  • PsiL50_H: Leaf water potential at 50% loss for high vegetation [MPa]
  • PsiL50_L: Leaf water potential at 50% loss for low vegetation [MPa]
  • PsiX50_H: Xylem water potential at 50% loss for high vegetation [MPa]
  • PsiX50_L: Xylem water potential at 50% loss for low vegetation [MPa]
  • Zs: Soil layer depths [mm]
  • row: Water density [kg/m³]

Returns

  • V: Water volume in each soil layer [mm]
  • O: Soil moisture in each soil layer [-]
  • OS: Surface soil moisture [-]
  • Lk: Leakage at bedrock [mm/h]
  • Psi_s_H: Soil water potential for high vegetation [MPa]
  • Psi_s_L: Soil water potential for low vegetation [MPa]
  • Exwat_H: Maximum extractable water for high vegetation [mm]
  • Exwat_L: Maximum extractable water for low vegetation [mm]
  • Rd: Saturation excess runoff [mm/dth]
  • TE_L: Updated low vegetation transpiration [kg/m²s]
  • TE_H: Updated high vegetation transpiration [kg/m²s]
  • E_soil: Updated soil evaporation [kg/m²s]
  • dV_dt: Change in water volume [mm/dth]
  • WBalance_soil: Water balance check [mm/dth]
  • Psi_soil: Soil water potential profile [mm]
  • Ko: Hydraulic conductivity profile [mm/h]
source
UrbanTethysChloris.Water.water_vegetationMethod
water_vegetation(
    Rain::FT,
    E_veg::FT,
    In_veg_tm1::FT,
    Sp_In::FT,
    LAI::FT,
    SAI::FT,
    row::FT,
    dth::FT
) where {FT<:AbstractFloat}

Calculate canopy rain interception and water balance for vegetation.

Arguments

  • Rain: precipitation [mm/dth]
  • E_veg: evaporation from vegetation [kg/m²s]
  • In_veg_tm1: Interception from previous time step [mm]
  • Sp_In: Specific water interception capacity [mm]
  • LAI: Leaf area index [-]
  • SAI: Stem area index [-]
  • row: density of water [kg/m³]
  • dth: calculation time step [h]

Returns

  • q_runon_veg: Runoff [mm/dth]
  • In_veg: Interception [mm]
  • dIn_veg_dt: Change in interception [mm/dth]
  • WBalance_In_veg: Water balance check [mm/dth]
source

Model Components

Forcing Inputs

UrbanTethysChloris.ModelComponents.ForcingInputs.AnthropogenicInputsType
AnthropogenicInputs{FT<:AbstractFloat}

Building and anthropogenic inputs affecting the urban environment.

Fields

  • Tb::Vector{FT}: Interior building temperature [K]
  • Qf_canyon::Vector{FT}: Anthropogenic heat input into the canyon air [W/m²]
  • Qf_roof::Vector{FT}: Anthropogenic heat input above the roof [W/m²]
  • Waterf_canyonVeg::Vector{FT}: Water applied on vegetated ground surface area [mm/time step]
  • Waterf_canyonBare::Vector{FT}: Water applied on bare ground surface area [mm/time step]
  • Waterf_roof::Vector{FT}: Water applied on roof surface area [mm/time step]
source
UrbanTethysChloris.ModelComponents.ForcingInputs.ForcingInputSetType
ForcingInputSet{FT<:AbstractFloat} <: AbstractForcingInputSet{FT}

Forcing inputs for the Urban Tethys-Chloris model.

Fields

  • datetime::Vector{DateTime}: Date and time of the forcing inputs.
  • anthropogenic::AnthropogenicInputs{FT}: Anthropogenic inputs.
  • hvacschedule::HVACSchedule{FT}: HVAC schedule.
  • meteorological::MeteorologicalInputs{FT}: Meteorological inputs.
  • sunposition::SunPositionInputs{FT}: Sun position inputs.
source
UrbanTethysChloris.ModelComponents.ForcingInputs.HVACScheduleType
HVACSchedule{FT<:AbstractFloat}

HVAC schedule parameters that specify heat and humidity sources from building equipment and occupants.

Fields

  • Hequip: Sensible heat from equipment [W/m² building ground area]
  • Hpeople: Sensible heat from people [W/m² building ground area]
  • LEequip: Latent heat from equipment [W/m² building ground area]
  • LEpeople: Latent heat from people [W/m² building ground area]
  • AirConRoomFraction: Fraction of air conditioned space [-]
source
UrbanTethysChloris.ModelComponents.ForcingInputs.MeteorologicalInputsType
MeteorologicalInputs{FT<:AbstractFloat}

Meteorological inputs and derived atmospheric properties.

Fields

  • LWR_in::Vector{FT}: Atmospheric longwave radiation [W/m² horizontal surface]
  • SAB1_in::Vector{FT}: Component of direct incoming shortwave radiation [W/m² horizontal surface]
  • SAB2_in::Vector{FT}: Component of direct incoming shortwave radiation [W/m² horizontal surface]
  • SAD1_in::Vector{FT}: Component of diffuse incoming shortwave radiation [W/m² horizontal surface]
  • SAD2_in::Vector{FT}: Component of diffuse incoming shortwave radiation [W/m² horizontal surface]
  • Tatm::Vector{FT}: Air Temperature at atmospheric reference level [K]
  • Uatm::Vector{FT}: Wind speed at atmospheric reference level [m/s]
  • Pre::Vector{FT}: Air pressure [Pa]
  • Rain::Vector{FT}: Precipitation [mm]
  • rel_hum::Vector{FT}: Relative humidity [-]
  • datetime::Vector{DateTime}: Timestamps for the data
  • esat_Tatm::Vector{FT}: Vapor pressure saturation at Tatm [Pa]
  • ea::Vector{FT}: Vapor pressure [Pa]
  • q_atm::Vector{FT}: Specific humidity of air at reference height [-]
  • qSat_atm::Vector{FT}: Saturation specific humidity [-]
  • SW_dir::Vector{FT}: Direct incoming shortwave radiation [W/m² horizontal surface]
  • SW_diff::Vector{FT}: Diffuse incoming shortwave radiation [W/m² horizontal surface]
  • Zatm::FT: Atmospheric reference height [m]
  • Catm_CO2::FT: Atmospheric CO2 concentration [ppm]
  • Catm_O2::FT: Intercellular Partial Pressure Oxygen [ppm]
  • SunDSM_MRT::FT: Mean radiant temperature from sun [K]
  • cp_atm::Vector{FT}: Specific heat of air [J/kg K]
  • rho_atm::Vector{FT}: Dry air density at atmosphere [kg/m³]
source
UrbanTethysChloris.ModelComponents.ForcingInputs.SunPositionInputsType
SunPositionInputs{FT<:AbstractFloat}

Solar position and timing parameters.

Fields

  • t_bef::FT: Time before solar noon [h]
  • t_aft::FT: Time after solar noon [h]
  • theta_Z::Vector{FT}: Solar zenith angle [rad]
  • theta_n::Vector{FT}: Difference between solar azimuth angle and canyon orientation [rad]
  • zeta_S::Vector{FT}: Solar azimuth angle [rad]
  • TimeOfMaxSolAlt::Vector{FT}: Time of maximum solar altitude [h]
source

Parameters

UrbanTethysChloris.ModelComponents.Parameters.HeightDependentVegetationParametersType
HeightDependentVegetationParameters{FT<:AbstractFloat} <: AbstractParameters{FT}

Parameters for the location-specific (tree, roof, ground) vegetation.

Fields

  • LAI::FT: Leaf area index [-]

  • SAI::FT: Stem area index [-].

  • hc::FT: Canopy height [m].

  • h_disp::FT: [-].

  • d_leaf::FT: [-].

  • CASE_ROOT::Int: Type of root profile [-].

  • ZR95::FT: Root depth, 95th percentile [mm].

  • ZR50::FT: Root depth, 50th percentile [mm].

  • ZRmax::FT: Root depth, maximum [mm]

  • Rrootl::FT: Root length index [m root m-2 PFT].

  • PsiL50::FT: Water potential at 50 % of leaf hydraulic conductivity [MPa].

  • PsiX50::FT: Water potential at 50 % of xylem hydraulic conductivity and limit for water extraction from soil [MPa].

  • FI::FT: Intrinsec quantum efficency [µmol CO2 µmol-1 photons].

  • Do::FT: Empirical coefficient that expresses the value of vapor pressure deficit at which f(∆e) = 0.5 [Pa].

  • a1::FT: Empirical parameter linking net assimilaton AnCto stomatal conductanceg{s,CO2}` [-].

  • go::FT: Minimum/cuticular stomatal conductance [mol CO2 m-2 leaf s-1].

  • CT::Int: Photosyntesis Typology for Plants, Photosynthetic pathway C3 or C4 [-]

  • DSE::FT: Activation Energy - Plant Dependent, Activation Energy in Photosynthesis for Rubisco Capacity [kJ/mol].

  • Ha::FT: Entropy factor - Plant Dependent, Activation energy [kJ / mol K].

  • gmes::FT: Mesophyll conductance, not used [mol CO2/s m2].

  • rjv::FT: Scaling factor between Jmax and V_{c,max} [µmol equivalent µmol-1 CO2].

  • Kopt::FT: Optical depth of direct beam per unit plant area (?) [-].

  • Knit::FT: Canopy nitrogen decay coefficient [-].

  • Vmax::FT: Maximum Rubisco capacity at 25◦C leaf scale [µmol CO2 m-2 s-1]

  • mSl::FT: [-].

  • e_rel::FT: Relative Efficiency of the photosynthesis apparatus due to Age/Day-length [-].

  • e_relN::FT: Relative efficiency of the photosynthesis apparatus due to N limitations [-].

  • Psi_sto_00::FT: Soil water potential at the beginning of stomatal closure [MPa].

  • Psi_sto_50::FT: Soil water potential at 50 % stomatal closure [MPa].

  • Sl::FT: Specific leaf area of biomass [m^2 /gC] [-].

  • SPARTREE::Int: Tree root distribution: 1 = Tree roots can access all water in the soil

(imp, bare, veg) equally; 2 = If the tree crown is smaller than the combined vegetated and bare fraction, then the trees only transpire from these fractions. Otherwise, they also transpire from the impervious ground fraction. [-].

source
UrbanTethysChloris.ModelComponents.Parameters.UrbanGeometryParametersType
UrbanGeometryParameters{FT<:AbstractFloat} <: AbstractParameters{FT}

Parameters for the UrbanGeometry model component.

Fields

  • Height_canyon::FT: Height of urban canyon [m].
  • Width_canyon::FT: Ground width of urban canyon [m].
  • Width_roof::FT: Roof width of urban canyon [m].
  • Height_tree::FT: Tree height [m].
  • Radius_tree::FT: Tree radius (=1/4 fg,tree *Wcan) [m].
  • Distance_tree::FT: Distance of wall to tree trunk [m].
  • Hcan_max::FT: Maximum height of roughness elements (buildings) [m].
  • Hcan_std::FT: Standard deviation of roughness elements (buildings) [m].
  • trees::Bool: Easy switch to include (=1) and exclude (=0) trees in the urban canyon.
  • ftree::FT: Tree fraction along canyon axis
  • hcanyon::FT: Normalized height of urban canyon [-].
  • wcanyon::FT: Normalized ground width of urban canyon [-].
  • wroof::FT: Normalized roof width of urban canyon [-]
  • htree::FT: Normalized tree height [-]
  • radius_tree::FT: Normalized tree radius [-].
  • distance_tree::FT: Normalized distance of wall to tree trunk [-].
  • ratio::FT: Height-to-width ratio [-].
  • wcanyon_norm::FT: Normalized canyon width overall [-].
  • wroof_norm::FT: Normalized roof width overall [-].
source
UrbanTethysChloris.ModelComponents.Parameters.VegetationParametersType
VegetationParameters{FT<:AbstractFloat} <: AbstractParameters{FT}

Container for vegetation parameters for different urban surface components.

Fields

  • roof::HeightDependentVegetationParameters{FT}: Vegetation parameters for roof vegetation
  • ground::HeightDependentVegetationParameters{FT}: Vegetation parameters for ground-level vegetation
  • tree::HeightDependentVegetationParameters{FT}: Vegetation parameters for trees
source
UrbanTethysChloris.ModelComponents.Parameters.ThermalPropertiesType
ThermalProperties{FT<:AbstractFloat} <: AbstractParameters{FT}

Container for vegetation parameters for different urban surface components.

Fields

  • roof::LocationSpecificThermalProperties{FT}: Roof thermal properties
  • ground::LocationSpecificThermalProperties{FT}: Ground thermal properties
  • wall::LocationSpecificThermalProperties{FT}: Wall thermal properties
source
UrbanTethysChloris.ModelComponents.Parameters.OpticalPropertiesType
OpticalProperties{FT<:AbstractFloat} <: AbstractParameters{FT}

Container for optical properties for different urban surface components.

Fields

  • roof::VegetatedOpticalProperties{FT}: Roof optical properties
  • ground::VegetatedOpticalProperties{FT}: Ground optical properties
  • wall::SimpleOpticalProperties{FT}: Wall optical properties
  • tree::SimpleOpticalProperties{FT}: Tree optical properties
source
UrbanTethysChloris.ModelComponents.Parameters.VegetatedOpticalPropertiesType
VegetatedOpticalProperties{FT<:AbstractFloat} <: AbstractParameters{FT}

Optical properties specific to vegetated surfaces (roof and ground)

Fields

  • aveg::FT: Vegetation surface albedo (-)
  • aimp::FT: Impervious surface albedo (-)
  • abare::FT: Bare surface albedo (-), only used for ground
  • eveg::FT: Vegetation surface emissivity (-)
  • eimp::FT: Impervious surface emissivity (-)
  • ebare::FT: Bare surface emissivity (-), only used for ground
source
UrbanTethysChloris.ModelComponents.Parameters.BuildingEnergyModelParametersType
BuildingEnergyModelParameters{FT<:AbstractFloat} <: AbstractParameters{FT}

Container for all building energy model parameters.

Fields

  • indoor_optical::IndoorOpticalProperties{FT}: Indoor surface optical properties
  • thermal::ThermalBuilding{FT}: Building thermal properties
  • windows::WindowParameters{FT}: Window parameters
  • hvac::HVACParameters{FT}: HVAC system parameters
source
UrbanTethysChloris.ModelComponents.Parameters.IndoorOpticalPropertiesType
IndoorOpticalProperties{FT<:AbstractFloat} <: AbstractParameters{FT}

Optical properties for indoor building surfaces.

Fields

  • abc::FT: Albedo ceiling (-)
  • abw::FT: Albedo wall (-)
  • abg::FT: Albedo ground (-)
  • abm::FT: Albedo internal mass (-)
  • ec::FT: Emissivity ceiling (-)
  • eg::FT: Emissivity ground (-)
  • ew::FT: Emissivity wall (-)
  • em::FT: Emissivity internal mass (-)
source
UrbanTethysChloris.ModelComponents.Parameters.ThermalBuildingType
ThermalBuilding{FT<:AbstractFloat} <: AbstractParameters{FT}

Parameters specifying thermal properties and dimensions of building internal mass.

Fields

  • IntMassOn::Int: Include building internal mass in calculation (0/1)
  • FloorHeight::FT: Average floor height in building (m)
  • dzFloor::FT: Average thickness of floors in building (m)
  • dzWall::FT: Average thickness of walls in building (m)
  • lan_ground_floor::FT: Building ground thermal conductivity (W/m K)
  • cv_ground_floor::FT: Building ground volumetric heat capacity (J/m³ K)
  • lan_floor_IntMass::FT: Internal mass floor thermal conductivity (W/m K)
  • cv_floor_IntMass::FT: Internal mass floor volumetric heat capacity (J/m³ K)
  • lan_wall_IntMass::FT: Internal mass wall thermal conductivity (W/m K)
  • cv_wall_IntMass::FT: Internal mass wall volumetric heat capacity (J/m³ K)
source
UrbanTethysChloris.ModelComponents.Parameters.WindowParametersType
WindowParameters{FT<:AbstractFloat} <: AbstractParameters{FT}

Parameters for building windows.

Fields

  • WindowsOn::Int: Include windows in simulation (0/1)
  • GlazingRatio::FT: Window-to-wall ratio (-)
  • Uvalue::FT: U-value of windows (W/m² K)
  • lan_windows::FT: Thermal conductivity of windows (W/m K)
  • cv_glass::FT: Volumetric heat capacity of glass (J/m³ K)
  • dztot::FT: Total thickness of glass layers (m)
  • SHGC::FT: Solar heat gain coefficient (-)
  • SolarTransmittance::FT: Solar radiation transmittance through windows (-)
  • SolarAbsorptivity::FT: Fraction of solar radiation absorbed by window (-)
  • SolarAlbedo::FT: Window albedo (-)
source
UrbanTethysChloris.ModelComponents.Parameters.HVACParametersType
HVACParameters{FT<:AbstractFloat} <: AbstractParameters{FT}

Parameters for HVAC system.

Fields

  • ACon::Int: Enable air conditioning (0/1)
  • Heatingon::Int: Enable heating (0/1)
  • TsetpointCooling::FT: Cooling setpoint temperature (K)
  • TsetpointHeating::FT: Heating setpoint temperature (K)
  • RHsetpointCooling::FT: Cooling setpoint relative humidity (%)
  • RHsetpointHeating::FT: Heating setpoint relative humidity (%)
  • ACH::FT: Air changes per hour (1/h)
  • COPAC::FT: Coefficient of performance for AC (-)
  • COPHeat::FT: Coefficient of performance for heating (-)
  • f_ACLatentToQ::FT: Fraction of latent heat removed by AC that is condensed (-)
source
UrbanTethysChloris.ModelComponents.Parameters.SoilParametersType
SoilParameters{FT<:AbstractFloat} <: AbstractParameters{FT}

Container for soil parameters for different urban surface components.

Fields

  • roof::VegetatedSoilParameters{FT}: Roof soil parameters
  • ground::VegetatedSoilParameters{FT}: Ground soil parameters
  • Sp_In_T::FT: Specific water retained by a tree (mm m^2 VEG area m^-2 plant area)
source
UrbanTethysChloris.ModelComponents.Parameters.VegetatedSoilParametersType
VegetatedSoilParameters{FT<:AbstractFloat} <: AbstractParameters{FT}

Soil parameters specific to vegetated surfaces (roof and ground)

Fields

  • Pcla::FT: Fraction of clay in the soil (-)
  • Psan::FT: Fraction of sand in the soil (-)
  • Porg::FT: Fraction of organic material in the soil (-)
  • In_max_imp::FT: Maxiumum interception capacity of impervious area (mm)
  • In_max_ground::FT: Maxiumum interception capacity of ground under roof vegetation (mm)
  • In_max_underveg::FT: Maxiumum interception capacity of vegetated ground area (mm)
  • In_max_bare::FT: Maxiumum interception capacity of bare ground area (mm)
  • Sp_In::FT: specific water retained by a vegetated surface (mm m^2 VEG area m^-2 plant area)
  • Kimp::FT: Hydraulic conductivity of impervious area (mm/h)
  • Kfc::FT: Conductivity at field capacity (mm/h)
  • Phy::FT: Suction at the residual/hygroscopic water content (kPa)
  • SPAR::Int: SOIL PARAMETER TYPE 1-VanGenuchten 2-Saxton-Rawls
  • Kbot::FT: Conductivity at the bedrock layer (mm/h)
source
UrbanTethysChloris.ModelComponents.Parameters.PersonParametersType
PersonParameters{FT<:AbstractFloat} <: AbstractParameters{FT}

Parameters for the person in the urban canyon for MRT calculations.

Fields

  • PositionPx::FT: Position within canyon [m]
  • PositionPz::FT: Height of centre of person [m]
  • PersonWidth::FT: Horizontal radius of ellipse describing person (=hip width / 2) [-]
  • PersonHeight::FT: Vertical radius of ellipse describing person (= height / 2) [-]
  • HeightWind::FT: Height for wind speed to calculate OTC [m]
source
UrbanTethysChloris.ModelComponents.Parameters.ParameterSetType
ParameterSet{FT<:AbstractFloat} <: AbstractParameterSet{FT}

Parameters for the Urban Tethys-Chloris model.

Fields

  • building_energy::BuildingEnergyModelParameters{FT}: Parameters for the building energy model.
  • person::PersonParameters{FT}: Parameters for the person.
  • soil::SoilParameters{FT}: Parameters for the soil.
  • surfacefractions::SurfaceFractions{FT}: Parameters for the surface fractions.
  • thermal::ThermalProperties{FT}: Parameters for the thermal properties.
  • optical::OpticalProperties{FT}: Parameters for the optical properties.
  • urbangeometry::UrbanGeometryParameters{FT}: Parameters for the urban geometry.
  • vegetation::VegetationParameters{FT}: Parameters for the vegetation.
source