Title: | Set of Tools to Compute Various Climate Indices |
---|---|
Description: | Set of tools to compute metrics and indices for climate analysis. The package provides functions to compute extreme indices, evaluate the agreement between models and combine theses models into an ensemble. Multi-model time series of climate indices can be computed either after averaging the 2-D fields from different models provided they share a common grid or by combining time series computed on the model native grid. Indices can be assigned weights and/or combined to construct new indices. |
Authors: | BSC-CNS [aut, cph], Victòria Agudetse [cre], Nuria Perez-Zanon [aut] , An-Chi Ho [ctb], Nicolau Manubens [ctb], Alasdair Hunter [aut], Louis-Philippe Caron [ctb], Eva Rifà [ctb] |
Maintainer: | Victòria Agudetse <[email protected]> |
License: | GPL-3 |
Version: | 0.3.3 |
Built: | 2024-11-22 06:26:49 UTC |
Source: | CRAN |
This function computes the mean and the percentage of agreement between anomalies.
AnoAgree(ano, membersdim, na.rm = TRUE, ncores = NULL)
AnoAgree(ano, membersdim, na.rm = TRUE, ncores = NULL)
ano |
A multidimensional array. |
membersdim |
The dimension in which models are stored. |
na.rm |
A logical indicating whether missing values should be removed. If
|
ncores |
The number of cores to be used when computing the agreement. |
An array of one dimension less than the ano
object, except for
one dimensional arrays or vectors, for which an array of dimension 1 called
'var' is returned.
# Example with random sample: a <- NULL for(i in 1:20) { a <- c(a, rnorm(6)) } dim(a) <- c(lat = 2, lon = 3, var = 4, mod = 5) agree <- AnoAgree(ano = a, membersdim = which(names(dim(a)) == 'mod'), na.rm = TRUE, ncores = NULL) print(agree) a <- rnorm(6) agree <- AnoAgree(ano = a, membersdim = 1, na.rm = TRUE, ncores = NULL) print(agree)
# Example with random sample: a <- NULL for(i in 1:20) { a <- c(a, rnorm(6)) } dim(a) <- c(lat = 2, lon = 3, var = 4, mod = 5) agree <- AnoAgree(ano = a, membersdim = which(names(dim(a)) == 'mod'), na.rm = TRUE, ncores = NULL) print(agree) a <- rnorm(6) agree <- AnoAgree(ano = a, membersdim = 1, na.rm = TRUE, ncores = NULL) print(agree)
This function splits an array into a list as required by PlotLayout function from package "s2dv" when parameter 'special_args' is used. The function ArrayToList allows to add names to the elements of the list in two different levels, the 'list' or the 'sublist'.
ArrayToList(data, dim, level = "list", names = NULL)
ArrayToList(data, dim, level = "list", names = NULL)
data |
A multidimensional array. |
dim |
A character string indicating the name of the dimension to split or an integer indicating the position of the dimension. |
level |
A string character 'list' or 'sublist' indicating if it should be a list or a sublist. By default it creates a list. |
names |
A vector of character strings to name the list (if it is a single string, it would be reused) or a single character string to name the elements in the sublist. |
A list of arrays of the length of the dimension set in parameter 'dim'.
data <- array(1:240, c(month = 12, member = 5, time = 4)) # Create a list: datalist <- ArrayToList(data, dim = 'month', level = 'list', names = month.name) class(datalist) class(datalist[[1]]) str(datalist) # Create a sublist: datalist <- ArrayToList(data, dim = 'month', level = 'sublist', names = 'dots') class(datalist) class(datalist[[1]]) class(datalist[[1]][[1]]) str(datalist)
data <- array(1:240, c(month = 12, member = 5, time = 4)) # Create a list: datalist <- ArrayToList(data, dim = 'month', level = 'list', names = month.name) class(datalist) class(datalist[[1]]) str(datalist) # Create a sublist: datalist <- ArrayToList(data, dim = 'month', level = 'sublist', names = 'dots') class(datalist) class(datalist[[1]]) class(datalist[[1]][[1]]) str(datalist)
This function computes the t90p, t10p, cdd or rx5day indices from n-dimensional arrays.
Climdex( data, metric, threshold = NULL, base.range = NULL, dates = NULL, timedim = NULL, calendar = NULL, ncores = NULL )
Climdex( data, metric, threshold = NULL, base.range = NULL, dates = NULL, timedim = NULL, calendar = NULL, ncores = NULL )
data |
A numeric n-dimensional array containing daily maximum or minimum temperature, wind speed or precipitation amount. |
metric |
The metric to be computed, either 't90p', 't10p', 'Wx', 'cdd' or 'rx5day'. |
threshold |
For the 't90p' and 't10p' metrics, an array of the 90th/10th
percentiles must be included. This parameter can be computed with the
|
base.range |
The years used for the reference period. If NULL (by default), all years are used. |
dates |
A vector of dates with a calendar attributes. If NULL (by default), the 'time' attributes of parameter 'data' are considered. |
timedim |
An integer number indicating the position of the time dimension
in the parameter |
calendar |
A character indicating the calendar type. |
ncores |
The number of cores to be used when computing the index. |
A list of length 2:
$result
, an array with the same dimensions as the input array,
except for the temporal dimension which is renamed to 'year', moved
to the first dimension position and reduce to annual resolution.
$years
, a vector of the corresponding years.
David Bronaugh for the Pacific Climate Impacts Consortium (2015). climdex.pcic: PCIC Implementation of Climdex Routines. R package version 1.1-6. http://CRAN.R-project.org/package=climdex.pcic
##Example synthetic data: data <- 1:(2 * 3 * 372 * 1) dim(data) <- c(lon = 2, lat = 3, time = 372, model = 1) time <- c(seq(ISOdate(1900, 1, 1), ISOdate(1900, 1, 31), "day"), seq(ISOdate(1901, 1, 1), ISOdate(1901, 1, 31), "day"), seq(ISOdate(1902, 1, 1), ISOdate(1902, 1, 31), "day"), seq(ISOdate(1903, 1, 1), ISOdate(1903, 1, 31), "day"), seq(ISOdate(1904, 1, 1), ISOdate(1904, 1, 31), "day"), seq(ISOdate(1905, 1, 1), ISOdate(1905, 1, 31), "day"), seq(ISOdate(1906, 1, 1), ISOdate(1906, 1, 31), "day"), seq(ISOdate(1907, 1, 1), ISOdate(1907, 1, 31), "day"), seq(ISOdate(1908, 1, 1), ISOdate(1908, 1, 31), "day"), seq(ISOdate(1909, 1, 1), ISOdate(1909, 1, 31), "day"), seq(ISOdate(1910, 1, 1), ISOdate(1910, 1, 31), "day"), seq(ISOdate(1911, 1, 1), ISOdate(1911, 1, 31), "day")) metadata <- list(time = list(standard_name = 'time', long_name = 'time', calendar = 'gregorian', units = 'days since 1970-01-01 00:00:00', prec = 'double', dim = list(list(name = 'time', unlim = FALSE)))) attr(time, "variables") <- metadata attr(data, 'Variables')$dat1$time <- time thres <- rep(10, 31 * 2 * 3) dim(thres) <- c(jdays = 31, lon = 2, lat = 3, model = 1) str(thres) clim <- Climdex(data, metric = "t90p", threshold = thres) str(clim)
##Example synthetic data: data <- 1:(2 * 3 * 372 * 1) dim(data) <- c(lon = 2, lat = 3, time = 372, model = 1) time <- c(seq(ISOdate(1900, 1, 1), ISOdate(1900, 1, 31), "day"), seq(ISOdate(1901, 1, 1), ISOdate(1901, 1, 31), "day"), seq(ISOdate(1902, 1, 1), ISOdate(1902, 1, 31), "day"), seq(ISOdate(1903, 1, 1), ISOdate(1903, 1, 31), "day"), seq(ISOdate(1904, 1, 1), ISOdate(1904, 1, 31), "day"), seq(ISOdate(1905, 1, 1), ISOdate(1905, 1, 31), "day"), seq(ISOdate(1906, 1, 1), ISOdate(1906, 1, 31), "day"), seq(ISOdate(1907, 1, 1), ISOdate(1907, 1, 31), "day"), seq(ISOdate(1908, 1, 1), ISOdate(1908, 1, 31), "day"), seq(ISOdate(1909, 1, 1), ISOdate(1909, 1, 31), "day"), seq(ISOdate(1910, 1, 1), ISOdate(1910, 1, 31), "day"), seq(ISOdate(1911, 1, 1), ISOdate(1911, 1, 31), "day")) metadata <- list(time = list(standard_name = 'time', long_name = 'time', calendar = 'gregorian', units = 'days since 1970-01-01 00:00:00', prec = 'double', dim = list(list(name = 'time', unlim = FALSE)))) attr(time, "variables") <- metadata attr(data, 'Variables')$dat1$time <- time thres <- rep(10, 31 * 2 * 3) dim(thres) <- c(jdays = 31, lon = 2, lat = 3, model = 1) str(thres) clim <- Climdex(data, metric = "t90p", threshold = thres) str(clim)
Function to combine climate indices for multiple models through addition, subtraction, division or averaging, optionally applying weights to each index.
CombineIndices(indices, weights = NULL, operation = "mean")
CombineIndices(indices, weights = NULL, operation = "mean")
indices |
List of n-dimensional arrays with equal dimensions to be combined. |
weights |
Vector of weights for the indices, whose length is the same as
the list of parameter |
operation |
The operation for combining the indices, either |
An array of the same dimensions as one of the elements in the
parameter indices
.
a <- matrix(rnorm(6), 2, 3) b <- matrix(rnorm(6), 2, 3) comb_ind <- CombineIndices(indices = list(a, b), weights = c(2, 1), operation = "add") print(comb_ind) a <- rnorm(24) dim(a) <- c(lon = 2, lat = 3, mod = 4) b <- rnorm(24) dim(b) <- c(lon = 2, lat = 3, mod = 4) comb_ind <- CombineIndices(indices = list(a, b), weights = c(2, 1), operation = "add") print(comb_ind)
a <- matrix(rnorm(6), 2, 3) b <- matrix(rnorm(6), 2, 3) comb_ind <- CombineIndices(indices = list(a, b), weights = c(2, 1), operation = "add") print(comb_ind) a <- rnorm(24) dim(a) <- c(lon = 2, lat = 3, mod = 4) b <- rnorm(24) dim(b) <- c(lon = 2, lat = 3, mod = 4) comb_ind <- CombineIndices(indices = list(a, b), weights = c(2, 1), operation = "add") print(comb_ind)
This function computes daily anomalies from a vector containing the daily time series.
DailyAno(data, jdays = NULL, dates = NULL, calendar = NULL, na.rm = TRUE)
DailyAno(data, jdays = NULL, dates = NULL, calendar = NULL, na.rm = TRUE)
data |
A vector of daily data. |
jdays |
A vector of the corresponding day of the year. This vector must
be the same length as parameter |
dates |
If |
calendar |
A character indicating the calendar type. |
na.rm |
A logical indicating whether missing values should be removed. If
|
A vector of daily anomalies of the same length as parameter
data
.
# Time series in a vector example: data <- 1:10 jdays <- c(rep(1, 5), rep(2, 5)) daily_anomaly <- DailyAno(data = data, jdays = jdays, na.rm = TRUE) print(daily_anomaly)
# Time series in a vector example: data <- 1:10 jdays <- c(rep(1, 5), rep(2, 5)) daily_anomaly <- DailyAno(data = data, jdays = jdays, na.rm = TRUE) print(daily_anomaly)
This function computes the diurnal temperature indicator, defined as the number of days where the diurnal temperature variation exceeds the vulnerability threshold (defined as the mean(tmax -tmin) + 5 from the reference period).
DTRIndicator( tmax, tmin, ref, by.seasons = TRUE, dates = NULL, timedim = NULL, calendar = NULL, ncores = NULL )
DTRIndicator( tmax, tmin, ref, by.seasons = TRUE, dates = NULL, timedim = NULL, calendar = NULL, ncores = NULL )
tmax |
A numeric multidimensional array containing daily maximum temperature. |
tmin |
A numeric multidimensional array containing daily minimum
temperature. This array must be the same dimensions as |
ref |
An output list from the |
by.seasons |
If TRUE (by default), the DTR is computed for each season (December-January-February, March-April-May, June-July-August and September-October-November) seperately. If FALSE is specified, the montly mean DTR is computed. |
dates |
A vector of dates with a calendar attributes. If NULL (by default), the 'time' attributes of parameter 'tmax' and 'tmin' are considered. |
timedim |
An integer number indicating the position of the time dimension
in the parameters |
calendar |
A character indicating the calendar type. |
ncores |
The number of cores to be used when computing the index. |
A list of length 3:
$dtr.ref
, an array with the same dimensions as the input
data
, but with the time dimension reduce from daily to monthly
or seasonal resolution depending on the selected resolution in
by.season
.
$year
, a vector of the corresponding years.
$season
, a vector of the seasons or months corresponding to the
resolution selected in by.season
.
##Exmaple with synthetic data: tmax <- 1 : (2 * 3 * 730 * 1) dim(tmax) <- c(lon = 2, lat = 3, time = 730, model = 1) tmin <- (1 : (2 * 3 * 730 * 1)) - 1 dim(tmin) <- c(lon = 2, lat = 3, time = 730, model = 1) time <- seq(as.POSIXct("1900-01-01 12:00:00", tz = "", format = "%Y-%d-%m %H:%M:%S"), as.POSIXct("1901-31-12 18:00:00", tz = "", format = "%Y-%d-%m %H:%M:%S"), "day") time <- as.POSIXct(time, tz = "CET") metadata <- list(time = list(standard_name = 'time', long_name = 'time', calendar = 'noleap', units = 'days since 1970-01-01 00:00:00', prec = 'double', dim = list(list(name ='time', unlim = FALSE)))) attr(time, "variables") <- metadata attr(tmax, 'Variables')$dat1$time <- time attr(tmax, 'Variables')$common[[2]]$dim[[3]]$len = length(time) attr(tmax, 'Variables')$common[[2]]$dim[[3]]$vals <- time attr(tmin, 'Variables')$dat1$time <- time attr(tmin, 'Variables')$common[[2]]$dim[[3]]$len = length(time) attr(tmin, 'Variables')$common[[2]]$dim[[3]]$vals <- time a <- DTRRef(tmax, tmin, by.seasons = FALSE, ncores = NULL) aa <- DTRIndicator(tmax, tmin, ref = a, by.seasons = FALSE, ncores = NULL) str(aa) dim(aa$indicator)
##Exmaple with synthetic data: tmax <- 1 : (2 * 3 * 730 * 1) dim(tmax) <- c(lon = 2, lat = 3, time = 730, model = 1) tmin <- (1 : (2 * 3 * 730 * 1)) - 1 dim(tmin) <- c(lon = 2, lat = 3, time = 730, model = 1) time <- seq(as.POSIXct("1900-01-01 12:00:00", tz = "", format = "%Y-%d-%m %H:%M:%S"), as.POSIXct("1901-31-12 18:00:00", tz = "", format = "%Y-%d-%m %H:%M:%S"), "day") time <- as.POSIXct(time, tz = "CET") metadata <- list(time = list(standard_name = 'time', long_name = 'time', calendar = 'noleap', units = 'days since 1970-01-01 00:00:00', prec = 'double', dim = list(list(name ='time', unlim = FALSE)))) attr(time, "variables") <- metadata attr(tmax, 'Variables')$dat1$time <- time attr(tmax, 'Variables')$common[[2]]$dim[[3]]$len = length(time) attr(tmax, 'Variables')$common[[2]]$dim[[3]]$vals <- time attr(tmin, 'Variables')$dat1$time <- time attr(tmin, 'Variables')$common[[2]]$dim[[3]]$len = length(time) attr(tmin, 'Variables')$common[[2]]$dim[[3]]$vals <- time a <- DTRRef(tmax, tmin, by.seasons = FALSE, ncores = NULL) aa <- DTRIndicator(tmax, tmin, ref = a, by.seasons = FALSE, ncores = NULL) str(aa) dim(aa$indicator)
This function computes the mean diurnal temperature range (tmax - tmin).
DTRRef( tmax, tmin, by.seasons = TRUE, dates = NULL, timedim = NULL, calendar = NULL, na.rm = TRUE, ncores = NULL )
DTRRef( tmax, tmin, by.seasons = TRUE, dates = NULL, timedim = NULL, calendar = NULL, na.rm = TRUE, ncores = NULL )
tmax |
A numeric multidimensional array containing daily maximum temperature. |
tmin |
A numeric multidimensional array containing daily minimum temperature. |
by.seasons |
If TRUE (by default), the DTR is computed for each season (December-January-February, March-April-May, June-July-August and September-October-November) seperately. If FALSE is specified, the montly mean DTR is computed. |
dates |
A vector of dates with a calendar attributes. If NULL (by default), the 'time' attributes of parameter 'tmax' and 'tmin' are considered. |
timedim |
An integer number indicating the position of the time dimension
in the parameters |
calendar |
A character indicating the calendar type. |
na.rm |
A logical indicating whether missing values should be removed. If
|
ncores |
The number of cores to be used when computing the index. |
The function returns a reordered array with 'time' dimension in the
first position in the dtr.ref
label.
A list of length 2:
$dtr.ref
, an array with the same dimensions as the input
data
, but with the time dimension reduce from daily to monthly
or seasonal resolution depending on the selected resolution in
by.season
.
$season
, a vector of the season or months corresponding to the
resolution selected in by.season
.
##Exmaple with synthetic data: tmax <- 1:(2 * 3 * 365 * 1) dim(tmax) <- c(lon = 2, lat = 3, time = 365, model = 1) tmin <- (1:(2 * 3 * 365 * 1))-1 dim(tmin) <- c(lon = 2, lat = 3, time = 365, model = 1) time <- seq.Date(as.Date("1900-01-01", format = "%Y-%d-%m"), as.Date("1900-31-12", format = "%Y-%d-%m"), 1) time <- as.POSIXct(time, tz = "CET") metadata <- list(time = list(standard_name = 'time', long_name = 'time', calendar = 'noleap', units = 'days since 1970-01-01 00:00:00', prec = 'double', dim = list(list(name ='time', unlim = FALSE)))) attr(time, "variables") <- metadata attr(tmax, 'Variables')$dat1$time <- time attr(tmax, 'Variables')$common[[2]]$dim[[3]]$len = length(time) attr(tmax, 'Variables')$common[[2]]$dim[[3]]$vals <- time attr(tmin, 'Variables')$dat1$time <- time attr(tmin, 'Variables')$common[[2]]$dim[[3]]$len = length(time) attr(tmin, 'Variables')$common[[2]]$dim[[3]]$vals <- time a <- DTRRef(tmax, tmin, by.seasons = FALSE, ncores = NULL) str(a) tmax <- 1:(2 * 3 * 365 * 1) dim(tmax) <- c(2, 3, 365) tmin <- (1:(2 * 3 * 365 * 1))-1 dim(tmin) <- c(2, 3, 365) a <- DTRRef(tmax, tmin, by.seasons = FALSE, dates = time, timedim = 3, ncores = NULL) str(a)
##Exmaple with synthetic data: tmax <- 1:(2 * 3 * 365 * 1) dim(tmax) <- c(lon = 2, lat = 3, time = 365, model = 1) tmin <- (1:(2 * 3 * 365 * 1))-1 dim(tmin) <- c(lon = 2, lat = 3, time = 365, model = 1) time <- seq.Date(as.Date("1900-01-01", format = "%Y-%d-%m"), as.Date("1900-31-12", format = "%Y-%d-%m"), 1) time <- as.POSIXct(time, tz = "CET") metadata <- list(time = list(standard_name = 'time', long_name = 'time', calendar = 'noleap', units = 'days since 1970-01-01 00:00:00', prec = 'double', dim = list(list(name ='time', unlim = FALSE)))) attr(time, "variables") <- metadata attr(tmax, 'Variables')$dat1$time <- time attr(tmax, 'Variables')$common[[2]]$dim[[3]]$len = length(time) attr(tmax, 'Variables')$common[[2]]$dim[[3]]$vals <- time attr(tmin, 'Variables')$dat1$time <- time attr(tmin, 'Variables')$common[[2]]$dim[[3]]$len = length(time) attr(tmin, 'Variables')$common[[2]]$dim[[3]]$vals <- time a <- DTRRef(tmax, tmin, by.seasons = FALSE, ncores = NULL) str(a) tmax <- 1:(2 * 3 * 365 * 1) dim(tmax) <- c(2, 3, 365) tmin <- (1:(2 * 3 * 365 * 1))-1 dim(tmin) <- c(2, 3, 365) a <- DTRRef(tmax, tmin, by.seasons = FALSE, dates = time, timedim = 3, ncores = NULL) str(a)
This function returns the number of spells of more than
min.length
days which exceed or are below the given threshold
from daily data.
Extremes( data, threshold, op = ">", min.length = 6, spells.can.span.years = TRUE, max.missing.days = 5, dates = NULL, timedim = NULL, calendar = NULL, ncores = NULL )
Extremes( data, threshold, op = ">", min.length = 6, spells.can.span.years = TRUE, max.missing.days = 5, dates = NULL, timedim = NULL, calendar = NULL, ncores = NULL )
data |
A n-dimensional array containing daily data. |
threshold |
A n-dimensional array with the threshold to be/not to be
reach, usually given by the a percentile computed with the |
op |
The operator to use to compare data to threshold. |
min.length |
The minimum spell length to be considered. |
spells.can.span.years |
Whether spells can span years. |
max.missing.days |
Maximum number of NA values per time period. |
dates |
A vector of dates with a calendar attributes. If NULL (by default), the 'time' attributes of parameter 'data' are considered. |
timedim |
An integer number indicating the position of the time dimension
in the parameter |
calendar |
A character indicating the calendar type. |
ncores |
The number of cores to be used when computing the extreme. |
This routine compares data to the thresholds using the given
operator, generating a series of TRUE or FALSE values; these values are then
filtered to remove any sequences of less than min.length
days of TRUE
values. It then computes the lengths of the remaining sequences of TRUE values
(spells) and sums their lengths. The spells.can.spa .years
option
controls whether spells must always terminate at the end of a period, or
whether they may continue until the criteria ceases to be met or the end of
the data is reached. The default for fclimdex is FALSE.
A list of length 2:
$output1
, an array with the same dimensions as the original
data
, except the time dimension which is reduced to annual
resolution given a timeseries of maximum spell lengths for each year.
$year
, a vector indicating the corresponding years.
##Example synthetic data: data <- 1:(2 * 3 * 310 * 1) dim(data) <- c(time = 310, lon = 2, lat = 3, model = 1) time <- as.POSIXct(paste(sort(rep(1902:1911, 31)), 1, 1:31, sep = "-"), tz = "CET") metadata <- list(time = list(standard_name = 'time', long_name = 'time', calendar = 'noleap', units = 'days since 1970-01-01 00:00:00', prec = 'double', dim = list(list(name = 'time', unlim = FALSE)))) attr(time, "variables") <- metadata attr(data, 'Variables')$dat1$time <- time threshold <- Threshold(data, dates = NULL, base.range = NULL, qtiles = 0.9, ncores = NULL) res <- Extremes(data, threshold = threshold, op = ">", min.length = 6, spells.can.span.years = TRUE, max.missing.days = 5, ncores = NULL) str(res)
##Example synthetic data: data <- 1:(2 * 3 * 310 * 1) dim(data) <- c(time = 310, lon = 2, lat = 3, model = 1) time <- as.POSIXct(paste(sort(rep(1902:1911, 31)), 1, 1:31, sep = "-"), tz = "CET") metadata <- list(time = list(standard_name = 'time', long_name = 'time', calendar = 'noleap', units = 'days since 1970-01-01 00:00:00', prec = 'double', dim = list(list(name = 'time', unlim = FALSE)))) attr(time, "variables") <- metadata attr(data, 'Variables')$dat1$time <- time threshold <- Threshold(data, dates = NULL, base.range = NULL, qtiles = 0.9, ncores = NULL) res <- Extremes(data, threshold = threshold, op = ">", min.length = 6, spells.can.span.years = TRUE, max.missing.days = 5, ncores = NULL) str(res)
This auxiliary function returns the index of position of a region of longitudes in a given vector of longitudes.
Lon2Index(lon, lonmin, lonmax)
Lon2Index(lon, lonmin, lonmax)
lon |
vector of longitudes values. |
lonmin |
a numeric value indicating the minimum longitude of the region (understand as the left marging of the region). |
lonmax |
a numeric value indicating the maximum longitude of the region (understand as the right mariging of the region). |
the index of positions of all values inside the region in the vector lon.
lon <- 1 : 360 pos <- Lon2Index(lon, lonmin = -20, lonmax = 20) lon[pos] pos <- Lon2Index(lon, lonmin = 340, lonmax = 20) lon[pos] lon <- -180 : 180 pos <- Lon2Index(lon, lonmin = -20, lonmax = 20) lon[pos] pos <- Lon2Index(lon, lonmin = 340, lonmax = 20) lon[pos]
lon <- 1 : 360 pos <- Lon2Index(lon, lonmin = -20, lonmax = 20) lon[pos] pos <- Lon2Index(lon, lonmin = 340, lonmax = 20) lon[pos] lon <- -180 : 180 pos <- Lon2Index(lon, lonmin = -20, lonmax = 20) lon[pos] pos <- Lon2Index(lon, lonmin = 340, lonmax = 20) lon[pos]
This function selects the daily data corresponding to the specified season.
SeasonSelect(data, season, dates = NULL, timedim = NULL, calendar = NULL)
SeasonSelect(data, season, dates = NULL, timedim = NULL, calendar = NULL)
data |
A numeric multidimensional array containing daily data. |
season |
A charcater string indicating the season by the three months initials in capitals: 'DJF' for winter (summer), 'MAM' spring (autumn), 'JJA' for summer (winter) or 'SON' for autumn (spring) in the northern (southern) hemisphere. |
dates |
A vector of dates with a calendar attributes. If NULL (by default), the 'time' attributes of parameter 'data' are considered. |
timedim |
An integer number indicating the position of the time dimension
in the parameter |
calendar |
A character indicating the calendar type. |
A list of length 2:
$data
, a vector or array containing the daily values for the
selected season, with the same dimensions as data
input but the
'time' dimension reduce to the number of days corresponding to the
selected season.
$dates
, a vector of dates reduce to the number of days
corresponding to the selected season.
## Example with synthetic data: data <- 1:(2 * 3 * (366 + 365) * 2) dim(data) <- c(lon = 2, lat = 3, time = 366 + 365, model = 2) time <- seq(ISOdate(1903,1,1), ISOdate(1904,12,31), "days") time <- as.POSIXct(time, tz = "CET") metadata <- list(time = list(standard_name = 'time', long_name = 'time', calendar = 'noleap', units = 'days since 1970-01-01 00:00:00', prec = 'double', dim = list(list(name ='time', unlim = FALSE)))) attr(time, "variables") <- metadata attr(data, 'Variables')$dat1$time <- time attr(data, 'Variables')$dat2$time <- time attr(data, 'Variables')$common[[2]]$dim[[3]]$len = length(time) attr(data, 'Variables')$common[[2]]$dim[[3]]$vals <- time a <- SeasonSelect(data = data, season = 'JJA') str(a)
## Example with synthetic data: data <- 1:(2 * 3 * (366 + 365) * 2) dim(data) <- c(lon = 2, lat = 3, time = 366 + 365, model = 2) time <- seq(ISOdate(1903,1,1), ISOdate(1904,12,31), "days") time <- as.POSIXct(time, tz = "CET") metadata <- list(time = list(standard_name = 'time', long_name = 'time', calendar = 'noleap', units = 'days since 1970-01-01 00:00:00', prec = 'double', dim = list(list(name ='time', unlim = FALSE)))) attr(time, "variables") <- metadata attr(data, 'Variables')$dat1$time <- time attr(data, 'Variables')$dat2$time <- time attr(data, 'Variables')$common[[2]]$dim[[3]]$len = length(time) attr(data, 'Variables')$common[[2]]$dim[[3]]$vals <- time a <- SeasonSelect(data = data, season = 'JJA') str(a)
Subset a spatial region from spatial data giving a vector with the maximum and minimum of latitudes and longitudes of the selected region.
SelBox(data, lon, lat, region, londim = "lon", latdim = "lat", mask = NULL)
SelBox(data, lon, lat, region, londim = "lon", latdim = "lat", mask = NULL)
data |
An array with minimum two dimensions of latitude and longitude. |
lon |
Numeric vector of longitude locations of the cell centers of the
grid of |
lat |
Numeric vector of latitude locations of the cell centers of the
grid of |
region |
A vector of length four indicating the minimum longitude, the maximum longitude, the minimum latitude and the maximum latitude. |
londim |
A character string indicating the name of the longitudinal dimension. The default value is 'lon'. |
latdim |
A character string indicating the name of the latitudinal dimension. The default value is 'lat'. |
mask |
A matrix with the same spatial dimensions of |
A list of length 4:
$data
, an array with the same dimensions as the input
data
array, but with spatial dimension reduced to the selected
region
.
$lat
, a vector with the new corresponding latitudes for the
selected region
.
$lon
, a vector with the new corresponding longitudes for the
selected region
.
$mask
, if parameter mask
is supplied, an array with
reduced length of the dimensions to the selected region
.
Otherwise, a NULL element is returned.
# Example with synthetic data: data <- 1:(20 * 3 * 2 * 4) dim(data) <- c(lon = 20, lat = 3, time = 2, model = 4) lon <- seq(2, 40, 2) lat <- c(1, 5, 10) a <- SelBox(data = data, lon = lon, lat = lat, region = c(2, 20, 1, 5), londim = "lon", latdim = "lat", mask = NULL)
# Example with synthetic data: data <- 1:(20 * 3 * 2 * 4) dim(data) <- c(lon = 20, lat = 3, time = 2, model = 4) lon <- seq(2, 40, 2) lat <- c(1, 5, 10) a <- SelBox(data = data, lon = lon, lat = lat, region = c(2, 20, 1, 5), londim = "lon", latdim = "lat", mask = NULL)
Shift the longitudes of a data array. Only reasonable for global longitude shifting. It is useful for map plotting or aligning datasets.
ShiftLon(data, lon, westB, lon_dim = "lon", ncores = NULL)
ShiftLon(data, lon, westB, lon_dim = "lon", ncores = NULL)
data |
A named multidimensional array with at least 'lon_dim' dimension. |
lon |
A numeric vector of longitudes. The values are expected to be monotonic increasing. |
westB |
A number indicating the west boundary of the new longitudes. |
lon_dim |
A character string indicating the name of the longitude dimension in 'data'. The default value is 'lon'. |
ncores |
An integer indicating the number of cores used for computation. The default value is NULL (use only one core). |
A list of 2:
data |
Array of the shifted data with the same dimensions as parameter 'data'. |
lon |
The monotonic increasing new longitudes with the same length as parameter 'lon' and start at 'westB'. |
data <- array(data = 1:50, dim = c(lon = 360, lat = 181)) lon <- array(data = 0:359, dim = c(lon = 360)) lat <- -90:90 ## lat does not change shifted <- ShiftLon(data = data, lon = lon, westB = -180, ncores = 1) ## Not run: s2dv::PlotEquiMap(var = data, lon = lon, lat = lat, filled.continents = FALSE) s2dv::PlotEquiMap(var = shifted$data, lon = shifted$lon, lat = lat, filled.continents = FALSE) ## End(Not run)
data <- array(data = 1:50, dim = c(lon = 360, lat = 181)) lon <- array(data = 0:359, dim = c(lon = 360)) lat <- -90:90 ## lat does not change shifted <- ShiftLon(data = data, lon = lon, westB = -180, ncores = 1) ## Not run: s2dv::PlotEquiMap(var = data, lon = lon, lat = lat, filled.continents = FALSE) s2dv::PlotEquiMap(var = shifted$data, lon = shifted$lon, lat = lat, filled.continents = FALSE) ## End(Not run)
This function allows to subset (i.e. slice, take a chunk of) an array, in a
similar way as done in the function take()
in the package plyr. There
are two main snprovements:
First, the input array can have dimension
names, either in names(dim(x))
or in the attribute 'dimensions'. If
both exist, names(dim(x))
is prioritized. The dimensions to subset
along can be specified via the parameter along
either with integer
indices or either by their name.
Second, there are additional ways to
adjust which dimensions are dropped in the resulting array: either to drop
all, to drop none, to drop only the ones that have been sliced or to drop
only the ones that have not been sliced.
Subset(x, along, indices, drop = FALSE)
Subset(x, along, indices, drop = FALSE)
x |
A named multidimensional array to be sliced. It can have dimension
names either in |
along |
A vector with references to the dimensions to take the subset from: either integers or dimension names. |
indices |
A list of indices to take from each dimension specified in 'along'. If a single dimension is specified in 'along', it can be directly provided as an integer or a vector. |
drop |
Whether to drop all the dimensions of length 1 in the resulting array, none, only those that are specified in 'along', or only those that are not specified in 'along'. The possible values are: 'all' or TRUE, 'none' or FALSE, 'selected', and 'non-selected'. The default value is FALSE. |
An array with similar dimensions as the x
input, but with
trimmed or dropped dimensions.
#Example synthetic data: # Dimension has name already data <- 1:(2 * 3 * 372 * 1) dim(data) <- c(time = 372, lon = 2, lat = 3, model = 1) data_subset <- Subset(data, c('time', 'model'), list(1:10, TRUE), drop = 'selected') dim(data_subset) # Use attributes 'dimensions' data <- array(1:(2 * 3 * 372 * 1), dim = c(2, 3, 372, 1)) attributes(data)[['dimensions']] <- c('lat', 'lon', 'time', 'model') data_subset <- Subset(data, c('lon', 'lat'), list(1, 1), drop = TRUE) dim(data_subset)
#Example synthetic data: # Dimension has name already data <- 1:(2 * 3 * 372 * 1) dim(data) <- c(time = 372, lon = 2, lat = 3, model = 1) data_subset <- Subset(data, c('time', 'model'), list(1:10, TRUE), drop = 'selected') dim(data_subset) # Use attributes 'dimensions' data <- array(1:(2 * 3 * 372 * 1), dim = c(2, 3, 372, 1)) attributes(data)[['dimensions']] <- c('lat', 'lon', 'time', 'model') data_subset <- Subset(data, c('lon', 'lat'), list(1, 1), drop = TRUE) dim(data_subset)
This function computes the threshold based on a quantile value for each day of the year of the daily data input.
Threshold( data, dates = NULL, calendar = NULL, base.range = NULL, qtiles = 0.9, ncores = NULL, na.rm = FALSE )
Threshold( data, dates = NULL, calendar = NULL, base.range = NULL, qtiles = 0.9, ncores = NULL, na.rm = FALSE )
data |
A numeric n-dimensional array containing daily data. |
dates |
A vector of dates with a calendar attributes. If NULL (by default), the 'time' attributes of parameter 'data' is considered. |
calendar |
A character indicating the calendar type. |
base.range |
The years used for computing the threshold. |
qtiles |
Numeric vector with values between 0 and 1 indicating the quantiles to be computed. |
ncores |
The number of cores to be used when computing the threshold. |
na.rm |
A logical value. If TRUE, any NA and NaN's are removed before the quantiles are computed (default as FALSE). |
An array with similar dimensions as the data
input, but without
'time' dimension, and a new 'jdays' dimension.
##Example synthetic data: data <- 1:(2 * 3 * 372 * 1) dim(data) <- c(time = 372, lon = 2, lat = 3, model = 1) time <- as.POSIXct(paste(sort(rep(1900:1911, 31)), 1, 1:31, sep = "-"), tz = "CET") metadata <- list(time = list(standard_name = 'time', long_name = 'time', calendar = 'noleap', units = 'days since 1970-01-01 00:00:00', prec = 'double', dim = list(list(name = 'time', unlim = FALSE)))) attr(time, "variables") <- metadata attr(data, 'Variables')$dat1$time <- time a <- Threshold(data, dates = NULL, base.range = NULL, qtiles = 0.9, ncores = NULL) str(a)
##Example synthetic data: data <- 1:(2 * 3 * 372 * 1) dim(data) <- c(time = 372, lon = 2, lat = 3, model = 1) time <- as.POSIXct(paste(sort(rep(1900:1911, 31)), 1, 1:31, sep = "-"), tz = "CET") metadata <- list(time = list(standard_name = 'time', long_name = 'time', calendar = 'noleap', units = 'days since 1970-01-01 00:00:00', prec = 'double', dim = list(list(name = 'time', unlim = FALSE)))) attr(time, "variables") <- metadata attr(data, 'Variables')$dat1$time <- time a <- Threshold(data, dates = NULL, base.range = NULL, qtiles = 0.9, ncores = NULL) str(a)
This function computes the duration of a heat/cold wave as the number of consecutive days for which the maximum/minimum temperature is exceeding/below a threshold over a minimum number of days in month or seasonal resolution.
WaveDuration( data, threshold, op = ">", spell.length = 6, by.seasons = TRUE, dates = NULL, calendar = NULL, ncores = NULL )
WaveDuration( data, threshold, op = ">", spell.length = 6, by.seasons = TRUE, dates = NULL, calendar = NULL, ncores = NULL )
data |
A numeric n-dimensional array containing daily maximum or minimum temperature |
threshold |
An array with the threshold to be/not to be reach, usually
given by the 90th/10th percentiles for heat/cold waves computed with the
|
op |
A character ">" (by default) or ">=" for heat waves and "<" or "<=" for cold waves indicating the operator must be used to compare data to threshold. |
spell.length |
A number indicating the number of consecutive days with extreme temperature to be considered heat or cold wave. |
by.seasons |
If TRUE (by default), the wave duration is computed for each season (DJF/MAM/JJA/SON) separately. If FALSE is specified, the monthly wave duration is computed. |
dates |
A vector of dates including calendar attributes. If NULL (by default), the 'time' attributes of parameter 'data' is used. |
calendar |
A character indicating the calendar type. |
ncores |
The number of cores to be used when computing the wave duration. |
A list of length 2:
$result
, an array with the same dimensions as the input
data
, but with the time dimension reduce from daily to monthly
or seasonal resolution depending on the selected resolution in
by.season
.
$years
, a vector of the years and season/months corresponding
to the resolution selected in by.season
and temporal length of
the input data
.
##Example synthetic data: data <- 1:(2 * 3 * 31 * 5) dim(data) <- c(lon = 2, lat = 3, time = 31, model = 5) time <- as.POSIXct(paste(paste(1900, 1, 1:31, sep = "-"), paste(12, 0, 0.0, sep = ":")), tz = "CET") metadata <- list(time = list(standard_name = 'time', long_name = 'time', calendar = 'standard', units = 'days since 1970-01-01 00:00:00', prec = 'double', dim = list(list(name ='time', unlim = FALSE)))) attr(time, "variables") <- metadata attr(data, 'Variables')$dat1$time <- time threshold <- rep(40, 31) a <- WaveDuration(data, threshold, op = ">", spell.length = 6, by.seasons = TRUE, ncores = NULL) str(a)
##Example synthetic data: data <- 1:(2 * 3 * 31 * 5) dim(data) <- c(lon = 2, lat = 3, time = 31, model = 5) time <- as.POSIXct(paste(paste(1900, 1, 1:31, sep = "-"), paste(12, 0, 0.0, sep = ":")), tz = "CET") metadata <- list(time = list(standard_name = 'time', long_name = 'time', calendar = 'standard', units = 'days since 1970-01-01 00:00:00', prec = 'double', dim = list(list(name ='time', unlim = FALSE)))) attr(time, "variables") <- metadata attr(data, 'Variables')$dat1$time <- time threshold <- rep(40, 31) a <- WaveDuration(data, threshold, op = ">", spell.length = 6, by.seasons = TRUE, ncores = NULL) str(a)
This function performs square-root of the cosine of the latitude weighting on the given array.
WeightedCells(data, lat, lat_dim = "lat", method = "cos", ncores = NULL)
WeightedCells(data, lat, lat_dim = "lat", method = "cos", ncores = NULL)
data |
A numeric array with named dimensions, representing the data to be applied the weights. It should have at least the latitude dimension and it can have more other dimensions. |
lat |
A numeric vector or array with one dimension containing the latitudes (in degrees). |
lat_dim |
A character string indicating the name of the latitudinal dimension. The default value is 'lat'. |
method |
A character string indicating the type of weighting applied: 'cos' (cosine of the latitude) or 'sqrtcos' (square-root of the cosine of the latitude). The default value is 'cos'. |
ncores |
An integer indicating the number of cores to use for parallel computation. The default value is NULL. |
An array containing the latitude weighted data with same dimensions as parameter 'data'.
exp <- array(rnorm(1:30), dim = c(lat = 3, lon = 5, sdate = 2)) lat <- c(10, 15, 20) res <- WeightedCells(data = exp, lat = lat)
exp <- array(rnorm(1:30), dim = c(lat = 3, lon = 5, sdate = 2)) lat <- c(10, 15, 20) res <- WeightedCells(data = exp, lat = lat)
This function computes a spatial area-weighted average of n-dimensional arrays being possible to select a region and to add a mask to be applied when computing the average.
WeightedMean( data, lon, lat, region = NULL, mask = NULL, londim = "lon", latdim = "lat", na.rm = TRUE, ncores = NULL )
WeightedMean( data, lon, lat, region = NULL, mask = NULL, londim = "lon", latdim = "lat", na.rm = TRUE, ncores = NULL )
data |
A numeric array with named dimensions, representing the data to be applied the weights. It should have at least the latitude dimension and it can have more other dimensions. |
lon |
A numeric vector of longitude locations of the cell centers of the
grid of |
lat |
A numeric vector of latitude locations of the cell centers of the
grid of |
region |
A vector of length four indicating the minimum longitude, the maximum longitude, the minimum latitude and the maximum latitude of the region to be averaged. |
mask |
A matrix with the same spatial dimensions of |
londim |
A character string indicating the name of the longitudinal dimension. The default value is 'lon'. |
latdim |
A character string indicating the name of the latitudinal dimension. The default value is 'lat'. |
na.rm |
A logical value indicating whether missing values should be stripped before the computation proceeds, by default it is set to TRUE. |
ncores |
An integer indicating the number of cores to use for parallel computation. The default value is NULL. |
An array, matrix or vector containig the area-weighted average with
the same dimensions as data
, except for the spatial longitude and
latitude dimensions, which disappear.
# Example 1: data <- 1:(2 * 3 * 4 * 5) dim(data) <- c(lon = 2, lat = 3, time = 4, model = 5) lat <- c(1, 10, 20) lon <- c(1, 10) a <- WeightedMean(data = data, lon = lon, lat = lat, region = NULL) mask <- c(0, 1, 0, 1, 0, 1) dim(mask) <- c(lon = 2, lat = 3) a <- WeightedMean(data = data, lon = lon, lat = lat, mask = mask) region <- c(1, 10, 1, 10) a <- WeightedMean(data = data, lon = lon, lat = lat, region = region, mask = mask) # Example 2: data <- 1:(2 * 3 * 4) dim(data) <- c(lon = 2, lat = 3, time = 4) lat <- c(1, 10, 20) lon <- c(1, 10) a <- WeightedMean(data = data, lon = lon, lat = lat)
# Example 1: data <- 1:(2 * 3 * 4 * 5) dim(data) <- c(lon = 2, lat = 3, time = 4, model = 5) lat <- c(1, 10, 20) lon <- c(1, 10) a <- WeightedMean(data = data, lon = lon, lat = lat, region = NULL) mask <- c(0, 1, 0, 1, 0, 1) dim(mask) <- c(lon = 2, lat = 3) a <- WeightedMean(data = data, lon = lon, lat = lat, mask = mask) region <- c(1, 10, 1, 10) a <- WeightedMean(data = data, lon = lon, lat = lat, region = region, mask = mask) # Example 2: data <- 1:(2 * 3 * 4) dim(data) <- c(lon = 2, lat = 3, time = 4) lat <- c(1, 10, 20) lon <- c(1, 10) a <- WeightedMean(data = data, lon = lon, lat = lat)