gridTemp: Grids station temperatures

Description Usage Arguments Value Author(s) See Also Examples

View source: R/gridTemp.R

Description

Grids station interval temperature values, so that they can be used as MESH inputs. The gridding method is performed by the hydroTSM function hydrokrige, using the IDW (inverse distance weighting) algorithm. The gridding uses a basin-scale lapse rate, which is determined by fitting a linear model to the difference between each site's temperature and that of the lowest site, and the difference in elevation relative to the lowest site. The procedure is the same as used in the the function basinLapseRates. In effect, all site temperatures are converted to have the same elevation before gridding. After gridding, each temperature is raised to its specified elevation using the same lapse rate. Where there are only 1 or 2 stations with available air temperatures, the hourly x monthly lapse rates returned by basinLapseRates are used.

Usage

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
gridTemp(
  temp = NULL,
  source_file_name = "unknown",
  shed_raster = NULL,
  site_elev = NULL,
  lapse_rates = NULL,
  IDW_file = NULL,
  tmin = 223.15,
  tmax = 313.15,
  missing_value = NA_real_,
  quiet = TRUE,
  progress_bar = TRUE
)

Arguments

temp

Required. A list containing 3 elements: 1. the header meta data, 2. the column meta data, and 3. the air temperature values (in K). These values are returned automatically by the MESHr command read_tb0.

source_file_name

Required. The name of the original .tb0 source file. Default value is unknown. The name of the source file is written to the r2c file header.

shed_raster

Required. A RasterBrick object describing the MESH basin. This can be created using the MESHr command read_r2c_shed with the parameter as_rasters = TRUE.

site_elev

Required. A data frame of station elevations. Note that the first column must contain the station names (which must be the same as in the air temperatures), and the second column must contain the elevation (in m).

lapse_rates

Optional. If there are 2 or fewer air temperatures in any interval, then the lapse rate cannot be calculated. In this case, if the historical lapse_rates are specified, then they will be used. If they are not specified, then this function will terminate with an error message. So, if you are confident that your dataset always has at least 3 stations with non- missing values of air temperatures, then you can omit this parameter. Note that the historical lapse rates must be a data frame of 12 columns (monthly) and 24 rows (hourly) values as returned by the function basinLapseRates.

IDW_file

Required. Output file which holds gridded air temperatures for all time steps.

tmin

Required. The minimum permitted air temperature of the gridded (and lapsed) air temperatures. All values exceeding tmin will be set to this value. The default is 223.15 K, or -50 C.

tmax

Required. The maximum permitted air temperature of the gridded (and lapsed) air temperatures. All values exceeding tmax will be set to this value. The default is 313.15 K, or 40 C.

missing_value

Required. Value to be used if all values in an interval are missing. Default is NA_real_.

quiet

Optional. If TRUE (the default) messages are suppressed. If FALSE, the time interval and messages from each gridding are listed.

progress_bar

Optional. If TRUE (the default), a progress bar is displayed showing the completed fraction of the temp.

Value

If unsuccessful, returns FALSE. If successful, returns TRUE and the gridded temperatures values are written to the IDW_file. Note that each interval's temperatures are written as they is gridded. This saves on memory, but can be quite slow. Note that the air temperatures in the file are in K.

Author(s)

Kevin Shook

See Also

gridPrecip basinLapseRates read_r2c_shed read_tb0

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
## Not run: 
hourly_temp_file <- "Red_Deer_all_hourly_temp_new.tb0"
temp <- read_tb0(hourly_temp_file, values_only = FALSE, timezone = "Etc/GMT+7", NAvalue = -0.1)
source_file_name <- hourly_temp_file
shedfile <- "RedDeer_MESH_drainage_database.r2c"
shed_raster <- read_r2c_shed(shedfile, as_rasters = TRUE, values_only = TRUE)
elev_file <- "site_elevations.csv"
site_elev <- read.csv(elev_file, header = TRUE, stringsAsFactors = FALSE)
lapse_rates_file <- "RedDeerLapseRates.csv"
lapse_rates <- read.csv(lapse_rates_file, header = TRUE, stringsAsFactors = FALSE, row.names = 1)
IDW_file <- "RedDeerTemp.idw"
gridPrecip(temp = temp, source_file_name = source_file_name, 
shed_raster = shed_raster, site_elev = site_elev,
lapse_rates = lapse_rates)
## End(Not run)

CentreForHydrology/MESHr documentation built on Jan. 11, 2021, 8:34 p.m.