METSIM: Meteorology Simulator

MetSim is a meteorological simulator and forcing disaggregator for hydrologic modeling and climate applications. Metsim is based on MtClim and the preprocessor from version 4 of the VIC hydrologic model.

MetSim consists of 3 main modules that govern the operation of 3 major aspects of its operation:

1. Management of dataset preprocessing and IO

The MetSim object provides high level support for setting up jobs and infrastructure for running simulation/disaggregation steps. It is the main interface through which the other modules are accessed.

2. Simulation of meteorological forcings

The base implementation of the meteorological simulator is based off of the algorithms described in [1]. This component has been designed to be flexible in allowing for alternative implementations which may be specified during the setup of the MetSim object. The default implementation allows for the daily simulation of:

  • Mean daily temperature
  • Snow water equivalent (SWE)
  • Incoming shortwave radiation
  • Cloud cover fraction
  • Potential evapotranspiration
  • Vapor pressure

3. Disaggregation of daily simulation values to sub-daily timesteps

Daily data from given input or simulated via the forcings generation component of MetSim can be disaggregated down to sub-daily values at intervals specified in minutes (provided they divide evenly into 24 hours). The operation of these algorithms is also described in [1].

This documentation is a work in progress. If you don’t find what you’re looking for here, check out MetSim’s Github page.

Getting Started

Installation

MetSim itself is a pure Python package, but its dependencies are not. You should ensure that you have all of the required dependencies:

Then, install MetSim with pip:

$ pip install metsim

To run the test suite after installing MetSim, install py.test (pip install pytest) and run py.test --verbose.

Finally, you can install MetSim directly from the source if you desire to:

$ git clone https://github.com/UW-Hydro/MetSim.git
$ cd MetSim
$ python setup.py install
$ py.test --verbose

Basic Usage

MetSim provides a simple command line interface which is primarily operated via configuration files. For more information about the options available to be set in the configuration files see the Configuration Specifications page.

Once installed, MetSim can be used from the command line via:

ms /path/to/configuration [-v] [-n #]

Bracketed flags are optional; -v activates verbose mode to print messages about the status of a run, and -n activates parallelism. The number given after the -n flag is the number of processes to run. A good rule of thumb is to use one less process than the number of processsors (or threads) that the machine you are running on has.

References

[1](1, 2) Bohn, T. J., B. Livneh, J. W. Oyler, S. W. Running, B. Nijssen, and D. P. Lettenmaier, 2013a: Global evaluation of MTCLIM and related algorithms for forcing of ecological and hydrological models, Agr. Forest. Meteorol., 176, 38-49, doi:10.1016/j.agrformet.2013.03.003.
[2]Bristow, K.L., and G.S. Campbell, 1984. On the relationship between incoming solar radiation and daily maximum and minimum temperature. Agricultural and Forest Meteorology, 31:159-166.
[3]Running, S.W., R.R. Nemani, and R.D. Hungerford, 1987. Extrapolation of synoptic meteorological data in mountainous terrain and its use for simulating forest evaporation and photosynthesis. Canadian Journal of Forest Research, 17:472-483.
[4]Glassy, J.M., and S.W. Running, 1994. Validating diurnal climatology of the MT-CLIM model across a climatic gradient in Oregon. Ecological Applications, 4(2):248-257.
[5]Kimball, J.S., S.W. Running, and R. Nemani, 1997. An improved method for estimating surface humidity from daily minimum temperature. Agricultural and Forest Meteorology, 85:87-98.
[6]Thornton, P.E., and S.W. Running, 1999. An improved algorithm for estimating incident daily solar radiation from measurements of temperature, humidity, and precipitation. Agricultural and Forest Meteorology, 93:211-228.

Sitemap

Examples

Basics

Provided in the source are several examples that can help you to get started using MetSim. They are located in the examples directory. We will look at the example_nc.conf file. Its contents are:

# This is an example of an input file for MetSim
[MetSim]

# Time step in minutes
time_step = 60

# Forcings begin here (year/month/day:hour)
start = 1950/1/1:0

# Forcings end at this date (year/month/day)
stop = 1950/1/31

# Input specification
forcing = ./tests/data/test.nc
domain  = ./tests/data/domain.nc
state = ./tests/data/state_nc.nc
in_fmt = netcdf
domain_fmt = netcdf
state_fmt = netcdf

# Output specification
out_fmt = netcdf
out_dir = ./results
out_state = ./results/state.nc
out_prefix = forcing
out_precision = f8

# How to disaggregate
method = mtclim

[forcing_vars]
Prec = prec
Tmax = t_max
Tmin = t_min

[state_vars]
prec = prec
t_max = t_max
t_min = t_min

[domain_vars]
lat = lat
lon = lon
mask = mask
elev = elev

This is a minimal configuration file for MetSim, and contains 3 sections. The first section, [MetSim] describes some basic settings such as the locations of data and parameters used in calculations. For a complete description of the input format see Configuration Specifications. The key things to note in this section are the options specified under the # Input specification and # Output specification comment headers. The forcing and domain options refer to the two types of required input, and the in_format and out_format options tell MetSim how they should be treated.

The second and third sections ([forcing_vars] and [state_vars]) describe the variables in the datasets provided in the forcing and state options of the first section. The left side of the assignment is the name of the variable given in the forcing dataset, while the right hand side is the name the variable should be given within MetSim. Note that the variables shown here are the minimum required set to run the forcing generation. The names given on the right hand side are also important to name correctly, as they are referenced internally. If you are unsure what variable names are used internally see the Configuration Specifications page for a full breakdown.

To run this example from the command line, once you have installed MetSim, use the following command:

ms path/to/example_nc.conf --verbose

This will run MetSim and disaggregate to hourly data, and write out the results in a NetCDF file located in the directory specified under out_dir in the configuration file (here ./results). The addition of the --verbose flag provides some information back to you as MetSim runs. In the absence of this flag MetSim will quietly run in the background until finished, or some error has occurred.

Generating daily values

Daily values can be output by specifying a time_step of 1440 in the configuration file, such as the one shown in the previous section. This will prevent MetSim’s disaggregation routines from being run, and the results written out will be daily values.

Translating formats of daily values

Warning

This section only applies to daily input and output.

This section can be useful if you are interested in converting VIC format binary or ASCII forcing inputs into NetCDF format.

To configure this behavior, several configuration options will have to be set. First, the time_step variable must be set to 1440 to enable daily output. Then, the forcing_fmt and out_fmt variables should be specified. The final option that should be set is out_vars. This can be set to include only the variables you have in your input file, if you wish to not generate any new data, or it can be set to include any of the variables generated by the simulator specified in the method option.

Input Specifications

There are 3 required input data sources that must be specified in the configuration file or dictionary. Note that it is possible for a single file to be specified for all three sources, provided that it has all of the required data. For examples of the data see the tests/data directory within the MetSim code.

Input forcings

Specified as forcing in the configuration file. This can either be the path to a NetCDF file, or the path to a directory containing ASCII or binary data (in the VIC4 format). The input forcing data is used to provide the base forcing variables. The required variable data is minimum daily temperature, maximum daily temperature, and daily precipitation.

The variable names can be mapped via the configuration file in the forcing_vars section. For more information about how to set up your configuration file see the Configuration Specifications page.

Domain file

Specified as domain in the configuration file. The domain file provides information about the domain MetSim is to be run over. It is required to be a NetCDF file. The domain requires two variables to be valid.

First, is the mask variable, which provides information about which grid cells are valid to run MetSim on. Values that specify grid cells which should be processed are specified via a positive, finite number (one or greater). Cells which MetSim should ignore can be given as 0 or NaN.

Second, is the elev variable. This provides elevation data used for calculation of solar geometry. It should be specified in meters, and only needs to be given at sites which are marked to be processed via the mask variable.

It is important to ensure that all valid locations in mask have data in elev. Failure to ensure this will result in errors during runtime.

State file

The state file provides information about the history of each of the grid cells to be processed. There are four required variables.

The first two are daily minimum and daily maximum temperatures for the 90 days preceeding the start date specified in the configuration file. They should be specified as t_min and t_max respectively. Similarly precipitation should be given as prec. These variables are used to generate seasonal averages which are used in the calculation of shortwave and longwave radiation.

The final required variable is the initial snow water equivalent (SWE) for each grid cell. It should be named swe in the file.

Output Specifications

Attention

The time coordinate in MetSim’s output is local to the location of each cell! This means that for a single time slice in the NetCDF file all locations along a parallel (same latitude) will have the same solar geometry at that time.

The output variables that are available are dependent on the time step being used. There are two cases:

Daily Output

When time_step is set to 1440 in the configuration file, daily values are generated. The following variables are available for output at a daily time step:

  • t_min : Minimum temperature (also a required input value) (C)
  • t_max : Maximum temperature (also a required input value) (C)
  • prec : Precipitation (also a required input value) (mm/day)
  • swe : Snow water equivalent (mm)
  • vapor_pressure : Vapor pressure (kPa)
  • shortwave : Shortwave radiation (W/m^2)
  • tskc : Cloud cover fraction
  • pet : Potential evapotranpiration (mm/day)

Sub-daily Output

  • temp : Temperature (C)
  • prec : Precipitation (mm/timestep)
  • shortwave : Shortwave radiation (W/m^2)
  • vapor_pressure : Vapor pressure (kPa)
  • air_pressure : Air pressure (kPa)
  • rel_humid : Relative humidity
  • spec_humid : Specific humidity
  • longwave : Longwave radiation (W/m^2)
  • tsck : Cloud cover fraction
  • wind : Wind speed (only if given as an input) (m/s)

Configuration Specifications

This page documents the various options and parameters that can be set in the configuration file.

MetSim Section

Required Variables

time_step :: int: The timestep to disaggregate in minutes. If given as 1440 (number of minutes in a day), no disaggregation will occur. This value must divide 1440 evenly.

start :: str: The time to start simulation given in the format yyyy/mm/dd or yyyy/mm/dd:hh.

stop :: str: The time to end simulation given in the format yyyy/mm/dd.

forcing :: path: The path to the input forcing file(s). See the section on __forcing_vars__ for more details.

domain :: path: The path to the input domain file. See the section on __domain_vars__ for more details.

state :: path: The path to the input state file.

out_dir :: path: The location to write output to. If this path doesn’t exist, it will be created.

out_state :: path/filename.nc: The location to write state file to.

forcing_fmt :: str: A string representing the type of input files specified in the forcing entry. Can be one of the following: ascii, binary, netcdf, or data.

state_fmt :: str: A string representing the type of state file specified in the state entry. Can be either netcdf or data.

domain_fmt :: str: A string representing the type of state file specified in the domain entry. Can be either netcdf or data.

out_fmt:: str: A string representing the type of output to write to out_dir. Can be either netcdf, data, or ascii.

method :: str: A string representing the simulation methods to use. The current implementation only supports mtclim.

Optional Variables

out_prefix :: str: The output file base name. Defaults to forcing.

out_precision :: str: Precision to use when writing output. Defaults to f8. Can be either f4 or f8.

annual :: bool: Whether to chunk up the timeseries into years for processing. This option is useful to set for when you are limited on memory. Each year of output is written as {out_prefix}_{year} when active.

iter_dims :: list: The dimensions of input data to iterate over to accumulate sites. Defaults to ['lat', 'lon'].

verbose :: bool: Whether to print output to stdout. Should be set using the -v flag for command line usage. This can be set for scripting purposes, if desired. Set to 1 to print output; defaults to 0.

sw_prec_thresh :: float: Minimum precipitation threshold to take into account when simulating incoming shortwave radiation. Defaults to 0.

mtclim_swe_corr :: bool: Whether to activate MtClim’s SWE correction algorithm. Default to False.

lw_cloud :: str: Type of cloud correction to longwave radiation to apply. Can be either DEFAULT or CLOUD_DEARDORFF. Defaults to CLOUD_DEARDORFF. Capitalization does not matter.

lw_type :: str: Type of longwave radiation parameterization to apply. Can be one of the following: DEFAULT, TVA, ANDERSON, BRUTSAERT, SATTERLUND, IDSO, or PRATA. Defaults to PRATA. Capitalization does not matter.

tdew_tol :: float: Convergence criteria for the iterative calculation of dewpoint temperature in MtClim. Defaults to 1e-6.

tmax_daylength_fraction :: float : Weight for calculation of time of maximum daily temperature. Must be between 0 and 1. Defaults to 0.67.

snow_crit_temp :: float: Critical temperature for snow to melt. Defaults to -6.0 C.

snow_melt_rate :: float: Melt rate when temperature is less than snow_crit_temp. Defaults to 0.042 cm/K.

rain_scalar :: float: Scale factor for calculation of cloudy sky transmittance. Defaults to 0.75, range should be between 0 and 1.

tday_coef :: float: Scale factor for calculation of daily mean temperature. Defaults to 0.45, range should be between 0 and 1.

lapse_rate :: float: Used to calculate atmospheric pressure. Defaults to 0.0065 K/m.

out_vars :: list : List of variables to write to output. Should be a list containing valid variables. The list of valid variables is dependent on which simulation method is used, as well as whether disaggregation is used. Defaults to ['temp', 'prec', 'shortwave', 'longwave', 'vapor_pressure', 'red_humid']. For more information about input and output variables see the Input Specifications page.

forcing_vars and state_vars section

The forcing_vars and state_vars sections are where you can specify which variables are in your input data, and the corresponding symbols which MetSim will recognize. The format of this section depends on the value given in the in_fmt entry in the MetSim section of the configuration file. See below for conventions for each input type.

netcdf and data

The in_vars section for NetCDF and xarray input acts as a mapping between the variable names in the input dataset to the variable names expected by MetSim. The format is given as netcdf_varname = metsim_varname. The minimum required variables given have metsim_varname``s corresponding to ``t_min, t_max, and prec; these variable names correspond to minimum daily temperature (Celcius), maximum daily temperature (Celcius), and precipitation (mm/day).

ascii

The in_vars section for ASCII input acts similarly to the NetCDF input format, except for one key point. Variables should be given as a tautology: the format is given as metsim_varname = metsim_varname. The order that the variables are given corresponds to the column numbers that they appear in the input files. The minimum required variables are t_min, t_max, and prec; these variable names correspond to minimum daily temperature (Celcius), maximum daily temperature (Celcius), and precipitation (mm/day).

binary

This section has an input style for binary files that mimics the VIC version 4 input style. Each line is specified as varname = scale cdatatype, where varname is the name that MetSim should use for the column, scale is a floating point scaling factor that should be applied after conversion from binary to floating point; the conversion applied by the scale is applied after the value in the input is converted from binary to the cdatatype specified for each variable. Valid cdatatype``s are ``signed and unsigned. signed values are interpreted as values which can be positive or negative, whereas unsigned values are interpreted as values that can only be greater than or equal to zero.

domain_vars section

The domain_vars section is where information about the domain file is given. Since the domain file is given as a NetCDF file this section has a similar format to that of the NetCDF input file format described above. That is, entries should be of the form netcdf_varname = metsim_varname. The minimum required variables have metsim_varname``s corresponding to ``lat, lon, mask, and elev; these variable names correspond to latitude, longitude, a mask of valid cells in the domain, and the elevation given in meters.

API reference

This page provides an auto-generated summary of metsim’s API. For more details and examples, refer to the relevant chapters in the main part of the documentation.

MetSim

The main object of the MetSim package. The MetSim object is used to set up and launch forcing generation and/or disaggregation routines.

The MetSim object uses a class dictionary to refer to the model setup, which can be modified after instantiation if necessary. Before calling run or launch on the instance it is required to call the load function to ensure that all of the required parameters have been set and that the input data is sufficient to provide the output specified.

class metsim.metsim.MetSim(params: dict)[source]

MetSim handles the distribution of jobs that write to a common file by launching muliple processes and queueing up their writeback so that work can be done while IO is happening.

launch()[source]

Farm out the jobs to separate processes

run()[source]

Kicks off the disaggregation and queues up data for IO

write(suffix='')[source]

Dispatch to the right function based on the configuration given

write_ascii(suffix)[source]

Write out as ASCII to the output file

write_netcdf(suffix: str)[source]

Write out as NetCDF to the output file

metsim.metsim.wrap_run(func: <built-in function callable>, loc: dict, params: dict, ds: xarray.core.dataset.Dataset, state: xarray.core.dataset.Dataset, disagg: bool, out_times: pandas.tseries.index.DatetimeIndex, year: str)[source]

Iterate over a chunk of the domain. This is wrapped so we can return a tuple of locs and df.

Parameters:
  • func (callable) – The function to call to do the work
  • loc (dict) – Some subset of the domain to do work on
  • params (dict) – Parameters from a MetSim object
  • ds (xr.Dataset) – Input forcings and domain
  • state (xr.Dataset) – State variables at the point of interest
  • disagg (bool) – Whether or not we should run a disagg routine
  • out_times (pd.DatetimeIndex) – Times to return (should be trimmed 1 day at each end from the given index)
  • year (str) – The year being run. This is used to add on extra times to make output smooth at endpoints if the run is chunked in time.
Returns:

A list of tuples arranged as (location, hourly_output, daily_output)

Return type:

results

Physics

physics

metsim.physics.atm_pres(elev: float, lr: float) → float[source]

Atmospheric pressure (Pa) as a function of elevation (m)

[1]Iribane, J.V., and W.L. Godson, 1981. Atmospheric Thermodynamics, 2nd Edition. D. Reidel Publishing Company, Dordrecht, The Netherlands (p. 168)
Parameters:
  • elev – Elevation in meters
  • lr – Lapse rate (K/m)
Returns:

Atmospheric pressure (Pa)

Return type:

pressure

metsim.physics.calc_pet(rad: <built-in function array>, ta: <built-in function array>, dayl: <built-in function array>, pa: float, dt: float = 0.2) → <built-in function array>[source]

Calculates the potential evapotranspiration for aridity corrections in calc_vpd(), according to Kimball et al., 1997

Parameters:
  • rad – daylight average incident shortwave radiation (W/m2)
  • ta – daylight average air temperature (deg C)
  • dayl – daylength (s)
  • pa – air pressure (Pa)
  • dt – offset for saturation vapor pressure calculation
Returns:

Potential evapotranspiration (cm/day)

Return type:

pet

metsim.physics.solar_geom[source]

Flat earth assumption

Parameters:
  • elev – Elevation in meters
  • lat – Latitude in decimal format
  • lr – Lapse rate in K/m
Returns:

(tiny_rad_fract, daylength, flat_potrad, tt_max0)

Return type:

sg

metsim.physics.svp[source]

Compute the saturated vapor pressure.

[2]Maidment, David R. Handbook of hydrology. McGraw-Hill Inc., 1992 Equation 4.2.2.
Parameters:
  • temp – Temperature (degrees Celsius)
  • a – (optional) parameter
  • b – (optional) parameter
  • c – (optional) parameter
Returns:

Saturated vapor pressure (Pa)

Return type:

svp

metsim.physics.svp_slope(temp: pandas.core.series.Series, a: float = 0.61078, b: float = 17.269, c: float = 237.3)[source]

Compute the gradient of the saturated vapor pressure as a function of temperature.

[3]Maidment, David R. Handbook of hydrology. McGraw-Hill Inc., 1992. Equation 4.2.3.
Parameters:temp – Temperature (degrees Celsius)
Returns:Gradient of d(svp)/dT.
Return type:dsvp_dT

MtClim

MTCLIM

metsim.methods.mtclim.calc_snowpack(df: pandas.core.frame.DataFrame, params: dict, snowpack: float = 0.0)[source]

Estimate snowpack as swe.

Parameters:
  • df – Dataframe with daily timeseries of precipitation and minimum temperature.
  • snowpack – (Optional - defaults to 0) Initial snowpack
metsim.methods.mtclim.calc_srad_hum(df: pandas.core.frame.DataFrame, sg: dict, elev: float, params: dict)[source]

Calculate shortwave, humidity

Parameters:
  • df – Dataframe containing daily timeseries
  • elev – Elevation in meters
  • params – A dictionary of parameters from the MetSim object
metsim.methods.mtclim.calc_t_air(df: pandas.core.frame.DataFrame, elev: float, params: dict)[source]

Calculate mean daily temperature

Parameters:
  • df – Dataframe with daily max and min temperatures
  • elev – Elevation in meters
  • params – Dictionary containing parameters from a MetSim object.
metsim.methods.mtclim.run(forcing: pandas.core.frame.DataFrame, params: dict, sg: dict, elev: float, swe: float)[source]

Run all of the mtclim forcing generation

Parameters:
  • forcing (pd.DataFrame) – The daily forcings given from input
  • params (dict) – Dictionary of parameters from a MetSim object
  • solar_geom (dict) – Solar geometry of the site
  • elev (float) – Elevation of the site being simulated
  • swe (float) – Initial snow water equivalent (SWE) for the site being simulated
Returns:

Dataframe of daily or subdaily forcings

Return type:

forcing

metsim.methods.mtclim.sw_hum_iter(df: pandas.core.frame.DataFrame, sg: dict, pa: float, pva: pandas.core.series.Series, dtr: pandas.core.series.Series, params: dict)[source]

Calculated updated values for dewpoint temperature and saturation vapor pressure.

Parameters:
  • df – Dataframe containing daily timeseries of cloud cover fraction, tfmax, swe, and shortwave radiation
  • sg – Solar geometry dictionary, calculated with metsim.physics.solar_geom.
  • pa – Air pressure in Pascals
  • pva – Vapor presure in Pascals
  • dtr – Daily temperature range
  • params – A dictionary of parameters from a MetSim object
Returns:

A tuple of dewpoint temperature and saturation vapor pressure

Return type:

(tdew, svp)

Disagg

Disaggregates daily data down to finer grained data using some heuristics

metsim.disaggregate.disaggregate(df_daily: pandas.core.frame.DataFrame, params: dict, solar_geom: dict, t_begin: list = None, t_end: list = None)[source]

Take a daily timeseries and scale it down to a finer time scale.

Parameters:
  • df_daily (pd.DataFrame) – Dataframe containing daily timeseries. Should be the result of one of the methods provided in the methods directory.
  • params (dict) – A dictionary containing the class parameters of the MetSim object.
  • solar_geom (dict) – A dictionary of solar geometry variables
  • t_begin (list) – List of t_min and t_max for day previous to the start of df_daily. None indicates no extension of the record.
  • t_end (list) – List of t_min and t_max for day after the end of df_daily. None indicates no extension of the record.
Returns:

A dataframe with sub-daily timeseries.

Return type:

df_disagg

metsim.disaggregate.longwave(air_temp: pandas.core.series.Series, vapor_pressure: pandas.core.series.Series, tskc: pandas.core.series.Series, params: dict)[source]

Calculate longwave. This calculation can be performed using a variety of parameterizations for both the clear sky and cloud covered emissivity. Options for choosing these parameterizations should be passed in via the params argument.

For more information about the options provided in this function see:

[1]Bohn, T.J., Livneh, B., Oyler, J.W., Running, S.W., Nijssen, B. and Lettenmaier, D.P., 2013. Global evaluation of MTCLIM and related algorithms for forcing of ecological and hydrological models. Agricultural and forest meteorology, 176, pp.38-49, doi:10.1016/j.agrformet.2013.03.003.
Parameters:
  • air_temp – Sub-daily temperature
  • vapor_pressure – Sub-daily vapor pressure
  • tskc – Daily cloud fraction
  • params – A dictionary of parameters, which contains information about which emissivity and cloud fraction methods to use.
Returns:

A sub-daily timeseries of the longwave radiation as well as a sub-daily timeseries of the cloud cover fraction.

Return type:

(lwrad, tskc)

metsim.disaggregate.prec(prec: pandas.core.series.Series, ts: float)[source]

Splits the daily precipitation evenly throughout the day Note: this returns only through to the beginning of the

last day. Final values are filled in using a forward fill in the top level disaggregate function
Parameters:
  • prec – Daily timeseries of precipitation
  • ts – Timestep to disaggregate down to
Returns:

A sub-daily timeseries of precipitation

Return type:

prec

metsim.disaggregate.pressure(temp: pandas.core.series.Series, elev: float, lr: float)[source]

Calculates air pressure.

Parameters:
  • temp – A sub-daily timeseries of temperature
  • elev – Elevation
  • lr – Lapse rate
Returns:

A sub-daily timeseries of air pressure (kPa)

Return type:

pressure

metsim.disaggregate.relative_humidity(vapor_pressure: pandas.core.series.Series, temp: pandas.core.series.Series)[source]

Calculate relative humidity from vapor pressure and temperature.

Parameters:
  • vapor_pressure – A sub-daily timeseries of vapor pressure
  • temp – A sub-daily timeseries of temperature
Returns:

A sub-daily timeseries of relative humidity

Return type:

rh

metsim.disaggregate.set_min_max_hour(disagg_rad: pandas.core.series.Series, n_days: int, ts: float, params: dict)[source]

Determine the time at which min and max temp is reached for each day.

Parameters:
  • disagg_rad – Shortwave radiation disaggregated to sub-daily timesteps.
  • n_days – The number of days being disaggregated
  • ts – Timestep used for disaggregation
  • params – A dictionary of class parameters of the MetSim object.
Returns:

A tuple containing 2 timeseries, corresponding to time of min and max temp, respectively

Return type:

(t_t_min, t_t_max)

metsim.disaggregate.shortwave(sw_rad: pandas.core.series.Series, daylength: pandas.core.series.Series, day_of_year: pandas.core.series.Series, tiny_rad_fract: <built-in function array>, params: dict)[source]

Disaggregate shortwave radiation down to a subdaily timeseries.

Parameters:
  • sw_rad – Daily incoming shortwave radiation
  • daylength – List of daylength time for each day of year
  • day_of_year – Timeseries of index of days since Jan-1
  • tiny_rad_fract – Fraction of the daily potential radiation during a radiation time step defined by SW_RAD_DT
  • params – Dictionary of parameters from the MetSim object
Returns:

A sub-daily timeseries of shortwave radiation.

Return type:

disaggrad

metsim.disaggregate.specific_humidity(vapor_pressure: pandas.core.series.Series, air_pressure: pandas.core.series.Series)[source]

Calculates specific humidity

Parameters:
  • vapor_pressure – A sub-daily timeseries of vapor pressure (kPa)
  • air_pressure – A sub-daily timeseries of air pressure (kPa)
Returns:

A sub-daily timeseries of specific humidity

Return type:

spec_humid

metsim.disaggregate.temp(df_daily: pandas.core.frame.DataFrame, df_disagg: pandas.core.frame.DataFrame, t_t_min: <built-in function array>, t_t_max: <built-in function array>, ts: float, t_begin: list = None, t_end: list = None)[source]

Disaggregate temperature using a Hermite polynomial interpolation scheme.

Parameters:
  • df_daily – A dataframe of daily values.
  • df_disagg – A dataframe of sub-daily values.
  • t_t_min – Times at which minimum daily temperatures are reached.
  • t_t_max – Times at which maximum daily temperatures are reached.
  • ts – Timestep for disaggregation
  • t_begin (list) – List of t_min and t_max for day previous to the start of df_daily. None indicates no extension of the record.
  • t_end (list) – List of t_min and t_max for day after the end of df_daily. None indicates no extension of the record.
Returns:

A sub-daily timeseries of temperature.

Return type:

temps

metsim.disaggregate.vapor_pressure(vp_daily: pandas.core.series.Series, temp: pandas.core.series.Series, t_t_min: <built-in function array>, n_out: int, ts: float)[source]

Calculate vapor pressure. First a linear interpolation of the daily values is calculated. Then this is compared to the saturated vapor pressure calculated using the disaggregated temperature. When the interpolated vapor pressure is greater than the calculated saturated vapor pressure, the interpolation is replaced with the saturation value.

Parameters:
  • vp_daily – Daily vapor pressure
  • temp – Sub-daily temperature
  • t_t_min – Timeseries of minimum daily temperature
  • n_out – Number of output observations
  • ts – Timestep to disaggregate down to
Returns:

A sub-daily timeseries of the vapor pressure

Return type:

vp

metsim.disaggregate.wind(wind: pandas.core.series.Series, ts: float)[source]

Wind is assumed constant throughout the day Note: this returns only through to the beginning of the

last day. Final values are filled in using a forward fill in the top level disaggregate function
Parameters:
  • wind – Daily timeseries of wind
  • ts – Timestep to disaggregate down to
Returns:

A sub-daily timeseries of wind

Return type:

wind