downscaleCV: Downscale climate data and reconstruct the temporal serie by...

View source: R/downscaleCV.R

downscaleCVR Documentation

Downscale climate data and reconstruct the temporal serie by splitting the data following a user-defined scheme

Description

Downscale climate data and reconstruct the temporal serie by splitting the data following a user-defined scheme. The statistical downscaling methods currently implemented are: analogs, generalized linear models (GLM) and Neural Networks (NN).

Usage

downscaleCV(
  x,
  y,
  method,
  sampling.strategy = "kfold.chronological",
  folds = 4,
  scaleGrid.args = NULL,
  simulate = FALSE,
  prepareData.args = list(global.vars = NULL, combined.only = TRUE, spatial.predictors
    = NULL, local.predictors = NULL, extended.predictors = NULL),
  condition = NULL,
  threshold = NULL,
  ...
)

Arguments

x

The input grid (admits both single and multigrid, see makeMultiGrid). It should be an object as returned by loadeR.

y

The observations dataset. It should be an object as returned by loadeR.

method

A string value. Type of transer function. Currently implemented options are "analogs", "GLM" and "NN".

sampling.strategy

Specifies a sampling strategy to define the training and test subsets. Possible values are "kfold.chronological" (the default), "kfold.random", "leave-one-year-out" and NULL. The sampling.strategy choices are next described:

  • "kfold.random" creates the number of folds indicated in the folds argument by randomly sampling the entries along the time dimension.

  • "kfold.chronological" is similar to "kfold.random", but the sampling is performed in ascending order along the time dimension.

  • "leave-one-year-out". This scheme performs a leave-one-year-out cross validation. It is equivalent to introduce in the argument folds a list of all years one by one.

  • NULL. The folds are specified by the user in the function parameter folds.

The first two choices will be controlled by the argument folds (see below)

folds

This arguments controls the number of folds, or how these folds are created (ignored if sampling.strategy = "leave-one-year-out"). If it is given as a fraction in the range (0-1), it splits the data in two subsets, one for training and one for testing, being the given value the fraction of the data used for training (i.e., 0.75 will split the data so that 75% of the instances are used for training, and the remaining 25% for testing). In case it is an integer value (the default, which is 4), it sets the number of folds in which the data will be split (e.g., folds = 10 for the classical 10-fold cross validation). Alternatively, this argument can be passed as a list, each element of the list being a vector of years to be included in each fold (See examples).

scaleGrid.args

A list of the parameters related to scale grids. This parameter calls the function scaleGrid. See the function definition for details on the parameters accepted.

simulate

A logic value indicating whether we want to simulate or not based on the GLM distributional parameters. Only relevant when perdicting with a GLM. Default to FALSE.

prepareData.args

A list with the arguments of the prepareData function. Please refer to prepareData help for more details about this parameter.

condition

Inequality operator to be applied considering the given threshold. "GT" = greater than the value of threshold, "GE" = greater or equal, "LT" = lower than, "LE" = lower or equal than. We only train with the days that satisfy the condition.

threshold

Numeric value. Threshold used as reference for the condition. Default is NULL. If a threshold value is supplied with no specification of the parameter condition. Then condition is set to "GE".

...

Optional parameters. These parameters are different depending on the method selected. Every parameter has a default value set in the atomic functions in case that no selection is wanted. Everything concerning these parameters is explained in the section Details of the function downscaleTrain. However, if wanted, the atomic functions can be seen here: analogs.train, glm.train and nn.train.

Details

The function relies on prepareData, prepareNewData, downscaleTrain, and downscalePredict. For more information please visit these functions. It is envisaged to allow for a flexible fine-tuning of the cross-validation scheme. It uses internally the transformeR helper dataSplit for flexible data folding. Note that the indices for data splitting are obtained using getYearsAsINDEX when needed (e.g. in leave-one-year-out cross validation), thus adequately handling potential inconsistencies in year selection when dealing with year-crossing seasons (e.g. DJF).

If the variable to downscale is the precipitation and it is a binary variable, then two temporal series will be returned:

  1. The temporal serie with binary values filtered by a threshold adjusted by the train dataset, see binaryGrid for more details.

  2. The temporal serie with the results obtained by the downscaling, without binary conversion process.

Missing data removal is recommended prior to multisite calibration.

According to the concept of cross-validation, a particular year should not appear in more than one fold (i.e., folds should constitute disjoint sets). For example, the choice fold =list(c(1988,1989), c(1989, 1990)) will raise an error, as 1989 appears in more than one fold.

Value

The reconstructed downscaled temporal serie.

Author(s)

J. Bano-Medina

See Also

downscaleTrain for training a downscaling model downscalePredict for prediction for a a test dataset with a trained model for downscaleR Wiki for downscaling seasonal forecasting and climate projections.

Other downscaling.functions: downscaleChunk(), downscalePredict(), downscaleTrain(), downscale()

Examples


require(transformeR)
require(climate4R.datasets)
data(NCEP_Iberia_hus850, NCEP_Iberia_ta850)
x <- makeMultiGrid(NCEP_Iberia_hus850, NCEP_Iberia_ta850)
x <- subsetGrid(x, years = 1985:1995)
# Loading predictands
data("VALUE_Iberia_pr")
y <- VALUE_Iberia_pr
y <- getTemporalIntersection(obs = y, prd = x, "obs")
x <- getTemporalIntersection(obs = y, prd = x, "prd")
# ... kfold in 3 parts equally divided ...
pred <- downscaleCV(x, y, folds = 3, sampling.strategy = "kfold.chronological",
                     scaleGrid.args = list(type = "standardize"),
                     method = "GLM", family = Gamma(link = "log"), condition = "GT", threshold = 0,
                     prepareData.args = list(
                     "spatial.predictors" = list(which.combine = getVarNames(x), v.exp = 0.9)))
# ... kfold by years ...
pred <- downscaleCV(x,y,sampling.strategy = "kfold.chronological",
                     method = "GLM", condition = "GT", threshold = 0,
                     scaleGrid.args = list(type = "standardize"),
                     folds = list(c(1985,1986,1987,1988),
                                  c(1989,1990,1991,1992),
                                  c(1993,1994,1995)))
# ... leave one year out  ...
pred <- downscaleCV(x,y,sampling.strategy = "leave-one-year-out",
                     method = "GLM", condition = "GT", threshold = 0,
                     scaleGrid.args = list(type = "standardize"))
# Reconstructing the downscaled serie in 3 folds with local predictors.
pred <- downscaleCV(x,y,folds = 3, sampling.strategy = "kfold.chronological",
                     scaleGrid.args = list(type = "standardize"),
                     method = "GLM", condition = "GT", threshold = 0,
                     prepareData.args = list("local.predictors" = list(vars = "hus@850", n = 4)))


SantanderMetGroup/downscaleR documentation built on July 4, 2023, 4:28 a.m.