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 daily 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
- 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]. The variables estimated are:
- Temperature
- Vapor pressure
- Relative and specific humidity
- Air pressure
- Cloud cover fraction
- Longwave radiation
- Shortwave radiation
- Precipitation
- Wind speed
For the “triangle” and “mix” methods of precipitation disaggregation, doumentation can be found here. This will eventually be superceded by a journal article that is currently in review [7].
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:
- Python 3.5 or greater
- xarray (0.10.9 or later)
- pandas (0.19.0 or later)
- numba (0.31.0 or later)
- netCDF4
- scipy
Then, install MetSim with pip or conda:
$ pip install metsim
or:
$ conda install -c conda-forge metsim
Alternatively, 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 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.
Warning
Users in environments where OpenMP is available may experience over-utilization of CPU resources, leading to lower performance. If you experience this issue try setting the OMP_NUM_THREADS environment variable to 1 before running MetSim.. This can be done in bash and similar shells by running export OMP_NUM_THREADS=1.
References¶
[1] | (1, 2) Bohn, T. J., B. Livneh, J. W. Oyler, S. W. Running, B. Nijssen, and D. P. Lettenmaier, 2013. Global evaluation of MTCLIM and related algorithms for forcing of ecological and hydrological models, Agricultural and Forest Meteorology, 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. |
[7] | Bohn, T. J., K. M. Whitney, G. Mascaro, and E. R. Vivoni, 2019. A deterministic approach for approximating the diurnal cycle of precipitation for large-scale hydrological simulations. Journal of Hydrometeorology, 20(2):297-317. doi: 10.1175/JHM-D-18-0203.1. |
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
# 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
# Output specification
out_dir = ./results
out_prefix = forcing
out_precision = f8
[chunks]
lat = 3
lon = 3
[forcing_vars]
prec = Prec
t_max = Tmax
t_min = Tmin
[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 the configuration page. 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 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 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 the following variables to be valid:
1. mask
: This 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
.
It is important to ensure that all valid locations in mask
have data in
elev
and any other variables. Failure to ensure this will result in
errors during runtime.
2. elev
: This provides elevation data (in m) used for calculation of solar
geometry. It only needs to be given at sites which are marked to be processed
via the mask
variable.
The next two variables are only needed if prec_type
= triangle
or
mix
in the input file:
3. dur
: This provides the climatological monthly storm event duration (in
minutes) used for disaggregating daily precipitation according to the
“triangle” method. Requires one value for each month (12).
4. t_pk
: This provides the climatological monthly time to storm peak (in
minutes starting from midnight) used for disaggregating daily precipitation to
sub-daily time scales using the “triangle” method. Requires one value for
each month (12).
For more information about the “triangle” method see this description. If you use this feature, please cite Bohn et al. (2019) as listed in the references.
A domain file for the CONUS+Mexico domain, at 0.0625 degree resolution, and
containing dur
and t_pk
values, is available here.
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.
Output Specifications¶
Attention
The time
coordinate in MetSim’s output is local to the location of each cell unless the utc_offset
is set to
True
! 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)vapor_pressure
: Vapor pressure (kPa)shortwave
: Shortwave radiation (W/m^2)tskc
: Cloud cover fractionpet
: Potential evapotranpiration (mm/day)wind
: Wind speed (only if given as an input) (m/s)
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 humidityspec_humid
: Specific humiditylongwave
: Longwave radiation (W/m^2)tsck
: Cloud cover fractionwind
: 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
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.
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
.
Optional Variables
output_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
.
time_grouper :: str
: Whether to chunk up the timeseries into pieces for
processing. This option is useful to set for when you are limited on
memory. Each chunk of output is written as {out_prefix}_{date_range}
when
active. Any valid pandas.TimeGrouper
string may be used (e.g. use ‘10AS’
for 10 year chunks).
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
.
rain_scalar :: float
: Scale factor for calculation of cloudy sky
transmittance. Defaults to 0.75
, range should be between 0
and
1
.
utc_offset :: bool
: Whether to use UTC timecode offsets for shifting
timeseries. Without this option all times should be considered local to
the gridcell being processed. Large domain runs probably want to set this
option to True
.
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
.
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']
.
prec_type :: str
: Type of precipitation disaggregation method to use. Can be
one of the following: uniform
, triangle
, or mix
. Defaults to
uniform
. Capitalization does not matter. Under uniform
method,
precipitation is disaggregated by dividing uniformly over all sub-daily
timesteps. Under triangle
the “triangle” method is employed whereby daily
precipitation is distributed assuming an isosceles triangle shape with peak and
width determined from two domain variables, t_pk
and dur
. Under
mix
, the “uniform” method is used on days when t_min
< 0 C, and
“triangle” is used on all other days; this hybrid method retains the improved
accuracy of “triangle” in terms of warm season runoff but avoids the biases
in snow accumulation that the “triangle” method sometimes yields due to fixed
event timing within the diurnal cycle of temperature. A domain file for the
CONUS+Mexico domain, containing the dur
and t_pk
parameters is
available at: <https://zenodo.org/record/1402223#.XEI-mM2IZPY>. For more
information about the “triangle” method see PtriangleMethod.pdf.
For more information about input and output variables see the Input Specifications page.
chunks section¶
The chunks
section describes how parallel computation should be grouped
in space. For example, to parallelize over 10 by 10 chunks of latitude and
longitude (with netcdf dimensions named lat
and lon
, respectively) you would use:
Alternatively, for an HRU based run chunked into 50 element jobs you would use:
As a general rule of thumb, try to evenly chunk the domain in such a way that the number of jobs to run is some multiple of the number of processors you wish to run on.
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 metsim_varname = netcdf_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 metsim_varname = netcdfvarname
. 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. If
prec_type
= triangle
or mix
, two additonal variables are required
including dur
and t_pk
for disaggregating daily precipitation according
to the “triangle” method.
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, domain_slice={})[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.
-
metsim.metsim.
add_prec_tri_vars
(domain)[source]¶ Check that variables for triangle precipitation method exist and have values that are within allowable ranges. Return these variables.
Parameters: domain – Dataset of domain variables for given location. Returns: - dur – Array of climatological monthly storm durations. [minutes]
- t_pk – Array of climatological monthly times to storm peaks. [minutes]
-
metsim.metsim.
chunk_domain
(chunks, dims)[source]¶ Return a dictionary of chunk slices that can be used to decompose a grid
-
metsim.metsim.
wrap_run_cell
(func: callable, params: dict, ds: xarray.core.dataset.Dataset, state: xarray.core.dataset.Dataset, disagg: bool, out_times: pandas.core.indexes.datetimes.DatetimeIndex)[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
- 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)
Returns: - df_complete (pd.DataFrame) – A dataframe with the disaggregated data in it
- df_base (pd.DataFrame) – A dataframe with the state data in it
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: numpy.array, ta: numpy.array, dayl: numpy.array, pa: float, dt: float = 0.2) → numpy.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
(temp: numpy.array, a: float = 0.61078, b: float = 17.269, c: float = 237.3)[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
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) → pandas.core.frame.DataFrame[source]¶ Take a daily timeseries and scale it down to a finer time scale.
Parameters: - df_daily – Dataframe containing daily timeseries. Should be the result of one of the methods provided in the methods directory.
- params – A dictionary containing the class parameters of the MetSim object.
- solar_geom – A dictionary of solar geometry variables
- t_begin – 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 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: numpy.array, vapor_pressure: numpy.array, tskc: numpy.array, params: dict) → numpy.array[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
Return type: lwrad
-
metsim.disaggregate.
prec
(prec: pandas.core.series.Series, t_min: pandas.core.series.Series, ts: float, params: dict, month_of_year: int)[source]¶ Distributes sub-daily precipitation either evenly (uniform) or with a triangular (triangle) distribution, depending upon the chosen method.
- Note: The uniform disaggregation 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. [mm]
- t_min – Daily timeseries of minimum daily temperature. [C]
- ts – Timestep length to disaggregate down to. [minutes]
- params – A dictionary of parameters, which contains information about which precipitation disaggregation method to use.
- month_of_year – Timeseries of index of month of year
Returns: A sub-daily timeseries of precipitation. [mm]
Return type: prec
-
metsim.disaggregate.
pressure
(temp: numpy.array, elev: float, lr: float) → numpy.array[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: numpy.array, temp: numpy.array) → numpy.array[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
(tiny_rad_fract: numpy.array, yday: numpy.array, n_days: int, ts: int, params: dict) → Tuple[numpy.array][source]¶ Determine the time at which min and max temp is reached for each day.
Parameters: - tiny_rad_fract – Array of fraction of shortwave radiation received at a shortened timestep. This should be calculated by metsim.physics.solar_geom.
- yday – Array of day of year for each simulated day.
- n_days – Number of days in run.
- ts – Timestep of run.
- params – Dictionary of parameters to use. Must contain utc_offset and tmax_daylength_fraction.
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: numpy.array, daylength: numpy.array, day_of_year: numpy.array, tiny_rad_fract: numpy.array, params: dict) → numpy.array[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: numpy.array, air_pressure: numpy.array) → numpy.array[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
(t_min: numpy.array, t_max: numpy.array, out_len: int, t_t_min: numpy.array, t_t_max: numpy.array, ts: int, t_begin: list = None, t_end: list = None) → numpy.array[source]¶ Disaggregate temperature using a Hermite polynomial interpolation scheme.
Parameters: - t_min – Timeseries of daily minimum temperatures.
- t_max – Timeseries of daily maximum temperatures.
- out_len – Length of the required output vector.
- 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 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 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.
tskc
(tskc: numpy.array, ts: int, params: dict) → numpy.array[source]¶ Disaggregate cloud fraction with uniform interpolation
Parameters: - tskc – Daily cloud fraction
- ts – Time step to disaggregate to (in minutes)
Returns: Sub-daily timeseries of cloud fraction
Return type: tskc
-
metsim.disaggregate.
vapor_pressure
(vp_daily: numpy.array, temp: numpy.array, t_t_min: numpy.array, n_out: int, ts: int) → numpy.array[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: numpy.array, ts: int, params: dict) → numpy.array[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 functionParameters: - wind – Daily timeseries of wind
- ts – Timestep to disaggregate down to
Returns: A sub-daily timeseries of wind
Return type: wind
What’s New¶
v.2.0.0¶
Enhancements¶
- Implemented UTC offsets, which puts all gridcell times in reference to UTC.
- Moved parallelism to dask, which allows for greater scalability and significantly less memory overhead.
Bug fixes¶
- Disallow timesteps > 6 hours, which raised errors.
- Raise error when t_min > t_max at beginning of runtime.
v.1.1.1¶
Enhancements¶
- Added option to disaggregate precipitation via a triangular hyetograph. (GH42). By Kristen Whitney and Theodore Bohn.
Bug fixes¶
- Fixed a bug where if iter_dims is not [lat, lon] the selected lat value that goes into solar_geom ends up as a list. The fix is also added for elevation and longitude, for redundancy. Fixes GH132. By Andrew Bennett.
v1.1.0¶
Enhancements¶
- Added option to use forcing start/stop dates to define run length (GH93). By Joe Hamman.
- Added option a flexible time grouper when chunking MetSim runs (GH93). By Joe Hamman.
- Improved configuration validation by checking for correctness of output variables (GH96) By Andrew Bennett.
- Added option to skip reading
swe
variable from state file if it is not going to be used by MtClim. (GH103). By Joe Hamman. - Added support for supplying a glob-like file path or multiple input forcing files (netCDF) (GH126). By Joe Hamman.
- Refactored
mtclim
anddisaggregate
functions to reduce interdependency and increase modularity. By Andrew Bennett. - Removed
swe
calculations. By Andrew Bennett.
Bug fixes¶
- Fixed bug where output files were not written with the appropriate calendar encoding attribute (GH97). By Joe Hamman.
- Fixed a bug where invalid timesteps were used in subdaily disaggregation. Added a clear error message explaining that subdaily timesteps must be evenly divisible into 24 hours and less than 6 hours in length. (GH110). By Joe Hamman.
- Fixed a bug during disaggregation when
t_min > t_max
. This now raises an exception. By Andrew Bennett.