API

Index

API documentation

PlantMeteo.ATMOSPHERE_COMPUTEDConstant
ATMOSPHERE_SELECT

List of variables that are by default removed from the table when using write_weather on a TimeStepTable{Atmosphere}.

source
PlantMeteo.AbstractAPIType
AbstractAPI

An abstract type for APIs. This is used to define the API to use for the weather forecast. You can get all available APIs using subtype(AbstractAPI).

source
PlantMeteo.AtmosphereType

Atmosphere structure to hold all values related to the meteorology / atmosphere.

Arguments

  • date<:AbstractDateTime = Dates.now(): the date of the record.
  • duration<:Period = Dates.Second(1.0): the duration of the time-step in Dates.Period.
  • T (°C): air temperature
  • Wind (m s-1): wind speed
  • Rh (0-1): relative humidity (can be computed using rh_from_vpd)
  • P = DEFAULTS.P (kPa): air pressure. The default value is at 1 atm, i.e. the mean sea-level

atmospheric pressure on Earth.

  • Precipitations = DEFAULTS.Precipitations (mm): precipitations from atmosphere (i.e. rain, snow, hail, etc.)
  • Cₐ = DEFAULTS.Cₐ (ppm): air CO₂ concentration
  • check = true: whether to check the validity of the input values.
  • e = vapor_pressure(T,Rh) (kPa): vapor pressure
  • eₛ = e_sat(T) (kPa): saturated vapor pressure
  • VPD = eₛ - e (kPa): vapor pressure deficit
  • ρ = air_density(T, P, constants.Rd, constants.K₀) (kg m-3): air density
  • λ = latent_heat_vaporization(T, constants.λ₀) (J kg-1): latent heat of vaporization
  • γ = psychrometer_constant(P, λ, constants.Cₚ, constants.ε) (kPa K−1): psychrometer "constant"
  • ε = atmosphere_emissivity(T,e,constants.K₀) (0-1): atmosphere emissivity
  • Δ = e_sat_slope(meteo.T) (0-1): slope of the saturation vapor pressure at air temperature
  • clearness::A = Inf (0-1): Sky clearness
  • Ri_SW_f::A = Inf (W m-2): Incoming short wave radiation flux
  • Ri_PAR_f::A = Inf (W m-2): Incoming PAR flux
  • Ri_NIR_f::A = Inf (W m-2): Incoming NIR flux
  • Ri_TIR_f::A = Inf (W m-2): Incoming TIR flux
  • Ri_custom_f::A = Inf (W m-2): Incoming radiation flux for a custom waveband

Notes

The structure can be built using only T, Rh, Wind and P. All other variables are optional and either let at their default value or automatically computed using the functions given in Arguments.

Examples

Atmosphere(T = 20.0, Wind = 1.0, P = 101.3, Rh = 0.65)
source
PlantMeteo.ConstantsType

Physical constants

The definition and default values are:

  • K₀ = -273.15: absolute zero (°C)
  • R = 8.314: universal gas constant ($J\ mol^{-1}\ K^{-1}$).
  • Rd = 287.0586: gas constant of dry air ($J\ kg^{-1}\ K^{-1}$).
  • Dₕ₀ = 21.5e-6: molecular diffusivity for heat at base temperature, applied in the integrated form of the Fick’s Law of diffusion ($m^2\ s^{-1}$). See eq. 3.10 from Monteith and Unsworth (2013).
  • Cₚ = 1013.0: Specific heat of air at constant pressure ($J\ K^{-1}\ kg^{-1}$), also known as efficiency of impaction of particles. See Allen et al. (1998), or Monteith and Unsworth (2013). NB: bigleaf R package uses 1004.834 intead.
  • ε = 0.622: ratio of molecular weights of water vapor and air. See Monteith and Unsworth (2013).
  • λ₀ = 2.501: latent heat of vaporization for water at 0 degree ($J\ kg^{-1}$).
  • σ = 5.670373e-08 Stefan-Boltzmann constant in ($W\ m^{-2}\ K^{-4}$).
  • Gbₕ_to_Gbₕ₂ₒ = 1.075: conversion coefficient from conductance to heat to conductance to water vapor.
  • Gsc_to_Gsw = 1.57: conversion coefficient from stomatal conductance to CO₂ to conductance to water vapor.
  • Gbc_to_Gbₕ = 1.32: conversion coefficient from boundary layer conductance to CO₂ to heat.
  • Mₕ₂ₒ = 18.0e-3 (kg mol-1): Molar mass for water.
  • J_to_umol = 4.57: Conversion factor from radiation in Joules m-2 s-1 (or W m-2) to μmol m-2 s-1.
  • PAR_fraction = 0.48: Fraction of shortwave radiation that is photosynthetically active radiation (PAR).

References

Allen, Richard G., Luis S. Pereira, Dirk Raes, et Martin J Fao Smith. 1998. « Crop evapotranspiration-Guidelines for computing crop water requirements-FAO Irrigation and drainage paper 56 » 300 (9): D05109.

Monteith, John, et Mike Unsworth. 2013. Principles of environmental physics: plants, animals, and the atmosphere. Academic Press.

source
PlantMeteo.OpenMeteoType
OpenMeteo()
OpenMeteo(
    vars=PlantMeteo.DEFAULT_OPENMETEO_HOURLY,
    forecast_server="https://api.open-meteo.com/v1/forecast",
    historical_server="https://archive-api.open-meteo.com/v1/era5",
    units=OpenMeteoUnits(),
    timezone="UTC",
    models=["auto"]
)

A type that defines the open-meteo.com API. No need of an API key as it is free. Please keep in mind that the API is distributed under the AGPL license, that it is not free for commercial use, and that you should use it responsibly.

Notes

The API wrapper provided by PlantMeteo is only working for the hourly data as daily data is missing some variables. The API wrapper is also not working for the "soil" variables as they are not consistant between forecast and historical data.

See also

to_daily

Arguments

  • vars: the variables needed, see here.
  • forecast_server: the server to use for the forecast, see

here. Default to https://api.open-meteo.com/v1/forecast.

  • historical_server: the server to use for the historical data, see

here. Default to https://archive-api.open-meteo.com/v1/era5.

  • start_archive::Dates.Day: the first day on which we have to get data from the historical archive instead of the forecast server,

data is at 25km resolution in the archive. Default to -150 days.

  • units::OpenMeteoUnits: the units used for the variables, see OpenMeteoUnits.
  • timezone: the timezone used for the data, see the list here.

Default to "UTC". This parameter is not checked, so be careful when using it.

  • models: the models to use for the forecast. Default to "["best_match"]". See OPENMETEO_MODELS for more details.

Troubleshooting

If you get an error when calling the API, try to decrease the value of start_archive to e.g.. 100 days before today (-Dates.Day(100)).

Details

The default variables are: "temperature2m", "relativehumidity2m", "precipitation", "surfacepressure", "windspeed10m", "shortwaveradiation", "directradiation", "diffuse_radiation".

Note that we don't download: "soiltemperature0cm", "soiltemperature6cm", "soiltemperature18cm", "soiltemperature54cm", "soilmoisture01cm", "soilmoisture13cm", "soilmoisture39cm", "soilmoisture927cm" and "soilmoisture27_81cm" by default as they are not consistant between forecast and hystorical data.

Sources

4.0 International (CC BY 4.0).

  • Copernicus Climate Change Service information 2022 (Hersbach et al., 2018).

Hersbach, H., Bell, B., Berrisford, P., Biavati, G., Horányi, A., Muñoz Sabater, J., Nicolas, J., Peubey, C., Radu, R., Rozum, I., Schepers, D., Simmons, A., Soci, C., Dee, D., Thépaut, J-N. (2018): ERA5 hourly data on single levels from 1959 to present. Copernicus Climate Change Service (C3S) Climate Data Store (CDS). (Updated daily), 10.24381/cds.adbb2d47

source
PlantMeteo.OpenMeteoUnitsType
OpenMeteoUnits(temperature_unit, windspeed_unit, precipitation_unit)
OpenMeteoUnits(;temperature_unit="celsius", windspeed_unit="ms", precipitation_unit="mm")

A type that defines the units used for the variabels when calling the open-meteo.com API.

Arguments

  • temperature_unit: the temperature unit, can be "celsius" or "fahrenheit". Default to "celsius".
  • windspeed_unit: the windspeed unit, can be "ms", "kmh", "mph", or "kn". Default to "ms".
  • precipitation_unit: the precipitation unit, can be "mm" or "inch". Default to "mm".

Examples

julia> units = OpenMeteoUnits("celsius", "ms", "mm")
OpenMeteoUnits("celsius", "ms", "mm")
source
PlantMeteo.TimeStepTableType
TimeStepTable(vars)

TimeStepTable stores variables values for each time step, e.g. weather variables. It implements the Tables.jl interface, so it can be used with any package that uses Tables.jl (like DataFrames.jl).

You can extend TimeStepTable to store your own variables by defining a new type for the storage of the variables. You can look at the Atmosphere type for an example implementation, or the Status type from PlantSimEngine.jl.

Examples

data = TimeStepTable(
    [
        Atmosphere(T = 20.0, Wind = 1.0, P = 101.3, Rh = 0.65),
        Atmosphere(T = 23.0, Wind = 1.5, P = 101.3, Rh = 0.60),
        Atmosphere(T = 25.0, Wind = 3.0, P = 101.3, Rh = 0.55)
    ]
)

# We can convert it into a DataFrame:
using DataFrames
df = DataFrame(data)

# We can also create a TimeStepTable from a DataFrame:
TimeStepTable(df)

# Note that by default it will use NamedTuple to store the variables
# for high performance. If you want to use a different type, you can
# specify it as a type parameter (if you want *e.g.* mutability or pre-computations):
TimeStepTable{Atmosphere}(df)
# Or if you use PlantSimEngine: TimeStepTable{Status}(df)
source
Base.setindex!Method
setindex!(row::TimeStepRow, nm::Symbol)
setindex!(row::TimeStepRow, i::Int)

Set the value of a variable in a TimeStepRow object.

source
PlantMeteo.WeatherMethod
Weather(data[, metadata])

Defines the weather, i.e. the local conditions of the Atmosphere for one or more time-steps. Each time-step is described using the Atmosphere structure, and the resulting structure is a TimeStepTable.

The simplest way to instantiate a Weather is to use a DataFrame as input.

The DataFrame should be formated such as each row is an observation for a given time-step and each column is a variable. The column names should match exactly the variables names of the Atmosphere structure:

See also

Examples

Example of weather data defined by hand (cumbersome):

w = Weather(
    [
        Atmosphere(T = 20.0, Wind = 1.0, P = 101.3, Rh = 0.65),
        Atmosphere(T = 23.0, Wind = 1.5, P = 101.3, Rh = 0.60),
        Atmosphere(T = 25.0, Wind = 3.0, P = 101.3, Rh = 0.55)
    ],
    (
        site = "Test site",
        important_metadata = "this is important and will be attached to our weather data"
    )
)

Weather is a TimeStepTable{Atmosphere}, so we can convert it into a DataFrame:

using DataFrames
df = DataFrame(w)

And then back into Weather to make a TimeStepTable{Atmosphere}:

Weather(df, (site = "My site",))

Of course it works with any DataFrame that has at least the required variables listed in Atmosphere.

source
PlantMeteo.add_transformations!Method
add_transformations!(df, trans, vars, fun; error_missing=false)

Add the fun transformations to the trans vector for the variables vars found in the df DataFrame.

Arguments

  • df: the DataFrame
  • trans: the vector of transformations (will be modified in-place)
  • vars: the variables to transform (can be a symbol, a vector of symbols or a pairs :var => :new_var)
  • fun: the function to apply to the variables
  • error_missing=true: if true, the function returns an error if the variable is not found. If false,

the variable is not added to trans.

source
PlantMeteo.air_densityMethod
air_density(Tₐ, P; check=true)
air_density(Tₐ, P, Rd, K₀; check=true)

ρ, the air density (kg m-3).

Arguments

  • Tₐ (Celsius degree): air temperature
  • P (kPa): air pressure
  • Rd (J kg-1 K-1): gas constant of dry air (see Foken p. 245, or R bigleaf package).
  • K₀ (Celsius degree): temperature in Celsius degree at 0 Kelvin
  • check (Bool): check if P is in the 85-110 kPa earth range

Note

Rd and K₀ are Taken from Constants if not provided.

References

Foken, T, 2008: Micrometeorology. Springer, Berlin, Germany.

source
PlantMeteo.atmosphere_emissivityMethod
atmosphere_emissivity(Tₐ,eₐ)

Emissivity of the atmoshpere at a given temperature and vapor pressure.

Arguments

  • Tₐ (°C): air temperature
  • eₐ (kPa): air vapor pressure
  • K₀ (°C): absolute zero

Examples

Tₐ = 20.0
VPD = 1.5
atmosphere_emissivity(Tₐ, vapor_pressure(Tₐ,VPD))

References

Leuning, R., F. M. Kelliher, DGG de Pury, et E.-D. SCHULZE. 1995. Leaf nitrogen, photosynthesis, conductance and transpiration: scaling from leaves to canopies ». Plant, Cell & Environment 18 (10): 1183‑1200.

source
PlantMeteo.check_day_completeMethod
check_day_complete(df)

Check that the weather table df has full days (24h) of data by summing their durations.

df must be a Tables.jl compatible table with a date and duration column. The date column must be a Dates.DateTime column, and the duration column must be a Dates.Period or Dates.CompoundPeriod column.

source
PlantMeteo.compute_dateFunction
compute_date(data, date_format, hour_format)

Compute the date column depending on several cases:

  • If it is already in data and is a DateTime, does nothing.
  • If it is a String, tries and parse it using a user-input DateFormat
  • If it is a Date, return it as is, or try to make it a DateTime if there's a column named

hour_start

Arguments

  • data: any Tables.jl compatible table, such as a DataFrame
  • date_format: a DateFormat to parse the date column if it is a String
  • hour_format: a DateFormat to parse the hour_start column if it is a String
source
PlantMeteo.compute_durationFunction
compute_duration(data, hour_format, duration)

Compute the duration column depending on several cases:

  • If it is already in the data, use the duration function to parse it into a Date.Period.
  • If it is not, but there's a column named hour_end and another one either called hour_start

or date, compute the duration from the period between hour_start (or date) and hour_end.

Arguments

  • data: any Tables.jl compatible table, such as a DataFrame
  • hour_format: a DateFormat to parse the hour_start and hour_end columns if they are Strings.
  • duration: a function to parse the duration column. Usually Dates.Day or Dates.Minute.
source
PlantMeteo.default_transformationMethod
default_transformation(df)

Return the default transformations to apply to the df DataFrame for the to_daily function. If the variable is not available, the transformation is not applied.

The default transformations are:

  • :date => (x -> unique(Dates.Date.(x))) => :date: the date is transformed into a Date object
  • :duration => sum => :duration: the duration is summed
  • :T => minimum => :Tmin: we use the minimum temperature for Tmin
  • :T => maximum => :Tmax: and the maximum temperature for Tmax
  • :T => mean => :T: and the average daily temperature (!= than average of Tmin and Tmax)
  • :Precipitations => sum => :Precipitations: the precipitations are summed
  • :Rh => mean => :Rh: the relative humidity is averaged
  • :Wind, :P, :Rh, :Cₐ, :e, :eₛ, :VPD, :ρ, :λ, :γ, :ε, :Δ, :clearness are

all averaged

  • [:Ri_SW_f, :duration] => ((x, y) -> sum(x .* Dates.toms.(y)) * 1.0e-9) => :Ri_SW_f: the irradiance is integrated, so its unit chaneges from W m⁻² to MJ m⁻² d⁻¹
  • All other irradiance variables are also integrated (see the code for details)
source
PlantMeteo.e_satMethod
e_sat(T)

Saturated water vapour pressure (es, in kPa) at given temperature T (°C). See Jones (1992) p. 110 for the equation.

source
PlantMeteo.fetch_openmeteoMethod
fetch_openmeteo(url, lat, lon, start_date, end_date, params::OpenMeteo)

Fetches the weather forecast from OpenMeteo.com and returns a tuple of:

  • a vector of Atmosphere
  • a NamedTuple of metadata (e.g. elevation, timezone, units...)
source
PlantMeteo.get_forecastMethod
get_forecast(params::OpenMeteo, lat, lon, period; verbose=true)

A function that returns the weather forecast from OpenMeteo.com

Arguments

  • lat: Latitude of the location
  • lon: Longitude of the location
  • period::Union{StepRange{Date, Day}, Vector{Dates.Date}}: Period of the forecast
  • verbose: If true, print more information in case of errors or warnings.
  • kwargs: Additional keyword arguments passed to get_forecast (not used in this method).

Examples

using PlantMeteo, Dates
lat = 48.8566
lon = 2.3522
period = [Dates.today(), Dates.today()+Dates.Day(3)]
params = OpenMeteo()
get_forecast(params, lat, lon, period)
source
PlantMeteo.get_weatherMethod
get_weather(lat, lon, period::Union{StepRange{Date, Day}, Vector{Dates.Date}}; api::DataType=OpenMeteo, sink=TimeStepTable, kwargs...)

Returns the weather forecast for a given location and time using a weather API.

Arguments

  • lat::Float64: Latitude of the location in degrees
  • lon::Float64: Longitude of the location in degrees
  • period::Union{StepRange{Date, Day}, Vector{Dates.Date}}: Period of the forecast
  • api::DataType=OpenMeteo: API to use for the forecast.
  • sink::DataType=TimeStepTable: Type of the output. Default is TimeStepTable, but it

can be any type that implements the Tables.jl interface, such as DataFrames.

  • kwargs...: Additional keyword arguments that are passed to the API

Details

We can get all available APIs using subtype(AbstractAPI). Please keep in mind that the default OpenMeteo API is not free for commercial use, and that you should use it responsibly.

Examples

using PlantMeteo, Dates
# Forecast for today and tomorrow:
period = [today(), today()+Dates.Day(1)]
w = get_weather(48.8566, 2.3522, period)
source
PlantMeteo.latent_heat_vaporizationMethod
latent_heat_vaporization(Tₐ,λ₀)
latent_heat_vaporization(Tₐ)

λ, the latent heat of vaporization for water (J kg-1).

Arguments

  • Tₐ (°C): air temperature
  • λ₀: latent heat of vaporization for water at 0 degree Celsius. Taken from Constants().λ₀

if not provided (see Constants).

source
PlantMeteo.new_namesMethod
new_names(args)

Get the new names of the columns after the transformation provided by the user.

source
PlantMeteo.next_valueFunction
next_value(row::TimeStepRow, var, next_index=1; default=nothing)

Return the value of var in the next row in the table, or default if there is no next row.

source
PlantMeteo.parse_hourFunction
parse_hour(h, hour_format=Dates.DateFormat("HH:MM:SS"))

Parse an hour that can be of several formats:

  • Time: return it as is
  • String: try to parse it using the user-input DateFormat
  • DateTime: transform it into a Time

Arguments

  • h: hour to parse
  • hour_format::DateFormat: user-input format to parse the hours

Examples

julia> using PlantMeteo, Dates;

As a string:

julia> PlantMeteo.parse_hour("12:00:00")
12:00:00

As a Time:

julia> PlantMeteo.parse_hour(Dates.Time(12, 0, 0))
12:00:00

As a DateTime:

julia> PlantMeteo.parse_hour(Dates.DateTime(2020, 1, 1, 12, 0, 0))
12:00:00
source
PlantMeteo.prepare_weatherMethod
prepare_weather(
    w;
    vars=setdiff(propertynames(w), ATMOSPHERE_COMPUTED),
    duration=Dates.Minute
)

Prepare the weather data for writing to a file. The function returns a DataFrame with the selected variables and the duration formated.

Arguments

  • w: a Tables.jl interfaced table, such as a TimeStepTable{Atmosphere} or a DataFrame
  • vars: a vector of variables to write (as symbols). By default, all variables are written except the ones that

can be recomputed (see ATMOSPHERE_COMPUTED). If nothing is given, all variables are written.

  • duration: the unit for formating the duration of the time steps. By default, it is Dates.Minute

Examples

using PlantMeteo, Dates

file = joinpath(dirname(dirname(pathof(PlantMeteo))),"test","data","meteo.csv")
w = read_weather(
    file,
    :temperature => :T,
    :relativeHumidity => (x -> x ./100) => :Rh,
    :wind => :Wind,
    :atmosphereCO2_ppm => :Cₐ,
    date_format = DateFormat("yyyy/mm/dd")
)

df = prepare_weather(w)
source
PlantMeteo.prev_rowFunction
prev_row(row::TimeStepRow, i=1)

Return the previous row in the table, or default if there is no previous row.

source
PlantMeteo.prev_valueFunction
prev_value(row::TimeStepRow, var, prev_index=1; default=nothing)

Return the value of var in the previous row in the table, or default if there is no previous row.

source
PlantMeteo.psychrometer_constantMethod
psychrometer_constant(P, λ, Cₚ, ε; check=true)
psychrometer_constant(P, λ; check=true)

γ, the psychrometer constant, also called psychrometric constant (kPa K−1). See Monteith and Unsworth (2013), p. 222.

Arguments

  • P (kPa): air pressure
  • λ ($J\ kg^{-1}$): latent heat of vaporization for water (see latent_heat_vaporization)
  • Cₚ (J kg-1 K-1): specific heat of air at constant pressure ($J\ K^{-1}\ kg^{-1}$)
  • ε (Celsius degree): temperature in Celsius degree at 0 Kelvin
  • check (Bool): check if P is in the 85-110 kPa earth range

Note

Cₚ, ε and λ₀ are taken from Constants if not provided.

Tₐ = 20.0

λ = latent_heat_vaporization(Tₐ, λ₀)
psychrometer_constant(100.0, λ)

References

Monteith, John L., et Mike H. Unsworth. 2013. « Chapter 13 - Steady-State Heat Balance: (i) Water Surfaces, Soil, and Vegetation ». In Principles of Environmental Physics (Fourth Edition), edited by John L. Monteith et Mike H. Unsworth, 217‑47. Boston: Academic Press.

source
PlantMeteo.read_weatherMethod
read_weather(
    file[,args...];
    date_format = DateFormat("yyyy-mm-ddTHH:MM:SS.s"),
    hour_format = DateFormat("HH:MM:SS"),
    duration=nothing
)

Read a meteo file. The meteo file is a CSV, and optionnaly with metadata in a header formatted as a commented YAML. The column names and units should match exactly the fields of Atmosphere, or the user should provide their transformation as arguments (args) with the DataFrames.jl form, i.e.:

  • :var_name => (x -> x .+ 1) => :new_name: the variable :var_name is transformed by the function x -> x .+ 1 and renamed to :new_name
  • :var_name => :new_name: the variable :var_name is renamed to :new_name
  • :var_name: the variable :var_name is kept as is

Note

The variables found in the file will be used as is if not transformed, and not recomputed from the other variables. Please check that all variables have the same units as in the Atmosphere structure.

Arguments

  • file::String: path to a meteo file
  • args...: A list of arguments to transform the table. See above to see the possible forms.
  • date_format = DateFormat("yyyy-mm-ddTHH:MM:SS.s"): the format for the DateTime columns
  • hour_format = DateFormat("HH:MM:SS"): the format for the Time columns (e.g. hour_start)
  • duration: a function to parse the duration column if present in the file. Usually Dates.Day or Dates.Minute.

If the column is absent, the duration will be computed using the hour_format and the hour_start and hour_end columns along with the date column.

Examples

using PlantMeteo, Dates

file = joinpath(dirname(dirname(pathof(PlantMeteo))),"test","data","meteo.csv")

meteo = read_weather(
    file,
    :temperature => :T,
    :relativeHumidity => (x -> x ./100) => :Rh,
    :wind => :Wind,
    :atmosphereCO2_ppm => :Cₐ,
    date_format = DateFormat("yyyy/mm/dd")
)
source
PlantMeteo.row_structMethod
row_struct(ts::TimeStepRow)

Get TimeStepRow in its raw format, e.g. the NamedTuple that stores the values, or the Atmosphere of values (or Status for PlantSimEngine.jl).

source
PlantMeteo.select_weatherFunction
select_weather(w, vars=setdiff(propertynames(w), ATMOSPHERE_COMPUTED))

Select the variables to write in the weather data. The function returns a DataFrame with the selected variables.

Arguments

  • w: a Tables.jl interfaced table, such as a TimeStepTable{Atmosphere} or a DataFrame
  • vars: a vector of variables to write (as symbols). By default, all variables are written except the ones that

can be recomputed (see ATMOSPHERE_COMPUTED). If nothing is given, all variables are written.

Examples

using PlantMeteo, Dates

file = joinpath(dirname(dirname(pathof(PlantMeteo))),"test","data","meteo.csv")
w = read_weather(
    file,
    :temperature => :T,
    :relativeHumidity => (x -> x ./100) => :Rh,
    :wind => :Wind,
    :atmosphereCO2_ppm => :Cₐ,
    date_format = DateFormat("yyyy/mm/dd")
)

df = select_weather(w)
source
PlantMeteo.standardize_columns!Method
standardize_columns!(::ToFileColumns, df)
standardize_columns!(::ToPlantMeteoColumns, df)

Standardize the column names of a DataFrame built upon Atmospheres to be compatible with standard file systems (CSV, databases...).

Arguments

  • df: a DataFrame built upon Atmospheres

Examples

using PlantMeteo, Dates, DataFrames

file = joinpath(dirname(dirname(pathof(PlantMeteo))),"test","data","meteo.csv")

df = read_weather(
    file,
    :temperature => :T,
    :relativeHumidity => (x -> x ./100) => :Rh,
    :wind => :Wind,
    :atmosphereCO2_ppm => :Cₐ,
    date_format = DateFormat("yyyy/mm/dd")
) |> DataFrame

df = standardize_columns!(ToFileColumns(), df)
source
PlantMeteo.sum_durationMethod
sum_duration(x::T) where {T<:AbstractVector{<:Dates.Period}}
sum_duration(x::T) where {T<:AbstractVector{<:Dates.CompoundPeriod}}

Sum the durations in x and return a Dates.Day object.

source
PlantMeteo.timesteps_durationsMethod
timesteps_durations(datetime::Dates.DateTime; verbose=true)

Duration in sensible units (e.g. 1 hour, or 1 day), computed as the duration between a step and the previous step. The first one is unknown, so we force it as the same as all (if unique), or the second one (if not) with a warning.

The function returns a Dates.CompoundPeriod because it helps finding a sensible default from a milliseconds period (e.g. 1 Hour or 1 Day).

Arguments

  • datetime::Vector{Dates.DateTime}: Vector of dates
  • verbose::Bool=true: If true, print a warning if the duration is not

constant between the time steps.

Examples

julia> timesteps_durations([Dates.DateTime(2019, 1, 1, 0), Dates.DateTime(2019, 1, 1, 1)])
2-element Vector{Dates.CompoundPeriod}:
 1 hour
 1 hour
source
PlantMeteo.to_dailyMethod
to_daily(df, args...)
to_daily(t::T, args...) where {T<:TimeStepTable{<:Atmosphere}}

Transform a DataFrame object or TimeStepTable{<:Atmosphere} with sub-daily time steps (e.g. 1h) to a daily time-step table.

Arguments

  • t: a TimeStepTable{<:Atmosphere} with sub-daily time steps (e.g. 1h)
  • args: a list of transformations to apply to the data, formates as for DataFrames.jl

Notes

Default transformations are applied to the data, and can be overriden by the user. The default transformations are:

  • :date => (x -> unique(Dates.Date.(x))) => :date: the date is transformed into a Date object
  • :duration => sum => :duration: the duration is summed
  • :T => minimum => :Tmin: we use the minimum temperature for Tmin
  • :T => maximum => :Tmax: and the maximum temperature for Tmax
  • :Precipitations => sum => :Precipitations: the precipitations are summed
  • :Rh => mean => :Rh: the relative humidity is averaged
  • :Wind, :P, :Rh, :Cₐ, :e, :eₛ, :VPD, :ρ, :λ, :γ, :ε, :Δ, :clearness are

all averaged

  • [:Ri_SW_f, :duration] => ((x, y) -> sum(x .* Dates.toms.(y)) * 1.0e-9) => :Ri_SW_f: the irradiance is integrated, so its unit chaneges from W m⁻² to MJ m⁻² d⁻¹
  • All other irradiance variables are also integrated (see the code for details)

Note that the default transformations can be overriden by the user, and that the default transformations are only applied if the variable is available.

Examples

using PlantMeteo, Dates
# Forecast for today and tomorrow:
period = [today(), today()+Dates.Day(1)]
w = get_weather(48.8566, 2.3522, period)
# Convert to daily:
w_daily = to_daily(w, :T => mean => :Tmean)
source
PlantMeteo.vapor_pressureMethod
vapor_pressure(Tₐ, Rh; check=true)

Vapor pressure (kPa) at given temperature (°C) and relative hunidity (0-1).

Arguments

  • Tₐ (Celsius degree): air temperature
  • Rh (0-1): relative humidity
  • check (Bool): if true, check that Rh is between 0 and 1

Examples

vapor_pressure(25.0, 0.4)
source
PlantMeteo.vpdMethod
vpd(Rh,eₛ)

Compute vapor pressure deficit (kPa) from the air relative humidity (0-1) and temperature (°C).

The computation simply uses vpd = eₛ - e.

Examples

vpd(0.4,25.0)
source
PlantMeteo.vpd_from_eMethod
vpd_from_e(e,Tₐ)

Compute vapor pressure deficit (kPa) from the air vapor pressure (kPa) and temperature (°C).

The computation simply uses vpd = eₛ - e.

Examples

vpd_from_e(1.5,25.0)
source
PlantMeteo.write_weatherMethod
write_weather(
    file, w; 
    vars=setdiff(propertynames(w), ATMOSPHERE_COMPUTED),
    duration=Dates.Minute
)

Write the weather data to a file.

Arguments

  • file: a String representing the path to the file to write
  • w: a TimeStepTable{Atmosphere}
  • vars: a vector of variables to write (as symbols). By default, all variables are written except the ones that

can be recomputed (see ATMOSPHERE_COMPUTED). If nothing is given, all variables are written.

  • duration: the unit for formating the duration of the time steps. By default, it is Dates.Minute.

Examples

using PlantMeteo, Dates

file = joinpath(dirname(dirname(pathof(PlantMeteo))),"test","data","meteo.csv")
w = read_weather(
    file,
    :temperature => :T,
    :relativeHumidity => (x -> x ./100) => :Rh,
    :wind => :Wind,
    :atmosphereCO2_ppm => :Cₐ,
    date_format = DateFormat("yyyy/mm/dd")
)

write_weather("meteo.csv", w)
source
Tables.getcolumnMethod
Tables.getcolumn(row::TimeStepRow, nm::Symbol)
Tables.getcolumn(row::TimeStepRow, nm::Int)

Get the value of a variable in a TimeStepRow object.

source