QUALYPSO: QUALYPSO

View source: R/QUALYPSO.r

QUALYPSOR Documentation

QUALYPSO

Description

Partition uncertainty in climate responses using an ANOVA applied to climate change responses.

Usage

QUALYPSO(Y, scenAvail, X = NULL, Xfut = NULL, iFut = NULL, listOption = NULL)

Arguments

Y

matrix nS x nY or array nG x nS x nY of climate projections.

scenAvail

data.frame nS x nEff with the nEff characteristics (e.g. type of GCM) for each of the nS scenarios. The number of characteristics nEff corresponds to the number of main effects that will be included in the ANOVA model.

X

(optional) predictors corresponding to the projections, e.g. time or global temperature. It can be a vector if the predictor is the same for all scenarios (e.g. X=2001:2100) or a matrix of the same size as Y if these predictors are different for the scenarios. By default, a vector 1:nY is created.

Xfut

(optional) nF values of the predictor over which the ANOVA will be applied. It must be a vector of values within the range of values of X. By default, it corresponds to X if X is a vector, 1:nY if X is NULL or a vector of 10 values equally spaced between the minimum and maximum values of X if X is a matrix.

iFut

index in 1:nF corresponding to a future predictor value . This index is necessary when Y is an array nG x nS x nY available for nG grid points. Indeed, in this case, we run QUALYPSO only for one future predictor. The first value defines the reference period or warming level.

listOption

(optional) list of options

  • args.smooth.spline: list of arguments to be passed to smooth.spline. The names attribute of args.smooth.spline gives the argument names (see do.call). The default option runs smooth.spline with spar=1.

  • typeChangeVariable: type of change variable: "abs" (absolute, value by default) or "rel" (relative).

  • ANOVAmethod: ANOVA method: "QUALYPSO" applies the method described in Evin et al. (2020), "lm" applies a simple linear model to estimate the main effects.

  • nBurn: if ANOVAmethod=="QUALYPSO", number of burn-in samples (default: 1000). If nBurn is too small, the convergence of MCMC chains might not be obtained.

  • nKeep: if ANOVAmethod=="QUALYPSO", number of kept samples (default: 2000). If nKeep is too small, MCMC samples might not represent correctly the posterior distributions of inferred parameters.

  • probCI: probability (in [0,1]) for the confidence intervals, probCI = 0.9 by default.

  • quantilePosterior: vector of probabilities (in [0,1]) for which we compute the quantiles from the posterior distributions quantilePosterior = c(0.005,0.025,0.05,0.1,0.25,0.33,0.5,0.66,0.75,0.9,0.95,0.975,0.995) by default.

  • climResponse: NULL by default. If it is provided, it must correspond to the outputs of fit.climate.response, i.e. a list with YStar [nS x nY], phiStar [nS x nF], etaStar [nS x nY], phi [nS x nF] and varInterVariability [scalar] if Y is a matrix [nS x nY], or a list with phiStar [nG x nS x nF], etaStar [nG x nS x nY], phi [nG x nS x nF] and varInterVariability vector of length nG if Y is an array [nG x nS x nY].

Value

List providing the results for each of the n values of Xfut if Y is a matrix or for each grid point if Y is an array, with the following fields:

  • CLIMATERESPONSE: list of climate change responses and corresponding internal variability. Contains phiStar (climate change responses), etaStar (deviation from the climate change responses as a result of internal variability), Ystar (change variable from the projections),and phi (fitted climate responses).

  • GRANDMEAN: List of estimates for the grand mean:

    • MEAN: vector of length n of means.

    • SD: vector of length n of standard dev. if ANOVAmethod=="QUALYPSO".

    • CI: matrix n x 2 of credible intervals of probability probCI given in listOption if ANOVAmethod=="QUALYPSO".

    • QUANT: matrix n x nQ of quantiles of probability quantilePosterior given in listOption if ANOVAmethod=="QUALYPSO".

  • MAINEFFECT: List of estimates for the main effects. For each main effect (GCM, RCM,..), each element of the list contains a list with:

    • MEAN: matrix n x nTypeEff

    • SD: matrix n x nTypeEff of standard dev. if ANOVAmethod=="QUALYPSO".

    • CI: array n x 2 x nTypeEff of credible intervals of probability probCI given in listOption if ANOVAmethod=="QUALYPSO".

    • QUANT: array n x nQ x nTypeEff of quantiles of probability quantilePosterior given in listOption if ANOVAmethod=="QUALYPSO".

  • CHANGEBYEFFECT: For each main effect, list of estimates for the mean change by main effect, i.e. mean change by scenario. For each main effect (GCM, RCM,..), each element of the list contains a list with:

    • MEAN: matrix n x nTypeEff

    • SD: matrix n x nTypeEff of standard dev. if ANOVAmethod=="QUALYPSO".

    • CI: array n x 2 x nTypeEff of credible intervals of probability probCI given in listOption if ANOVAmethod=="QUALYPSO".

    • QUANT: array n x nQ x nTypeEff of quantiles of probability quantilePosterior given in listOption if ANOVAmethod=="QUALYPSO".

  • EFFECTVAR: Matrix n x nTypeEff giving, for each time variability related to the main effects (i.e. variability between the different RCMs, GCMs,..).

  • CONTRIB_EACH_EFFECT: Contribution of each individual effect to its component (percentage), e.g. what is the contribution of GCM1 to the variability related to GCMs. For each main effect (GCM, RCM,..), each element of the list contains a matrix n x nTypeEff

  • RESIDUALVAR: List of estimates for the variance of the residual errors:

    • MEAN: vector of length n.

    • SD: vector of length n of standard dev. if ANOVAmethod=="QUALYPSO".

    • CI: matrix n x 2 of credible intervals of probability probCI given in listOption if ANOVAmethod=="QUALYPSO".

    • QUANT: matrix n x nQ of quantiles of probability quantilePosterior given in listOption if ANOVAmethod=="QUALYPSO".

  • INTERNALVAR: Internal variability (constant over time)

  • TOTALVAR: total variability, i.e. the sum of internal variability, residual variability and variability related to the main effects

  • DECOMPVAR: Decomposition of the total variability for each component

  • RESERR: differences between the climate change responses and the additive anova formula (grand mean + main effects)

  • Xmat: matrix of predictors

  • Xfut: future predictor values

  • paralType: type of parallelisation (Time or Grid)

  • namesEff: names of the main effects

  • Y: matrix of available combinations given as inputs

  • listOption: list of options used to obtained these results (obtained from QUALYPSO.check.option)

  • listScenarioInput: list of scenario characteristics (obtained from QUALYPSO.process.scenario)

Author(s)

Guillaume Evin

References

Evin, G., B. Hingray, J. Blanchet, N. Eckert, S. Morin, and D. Verfaillie (2020) Partitioning Uncertainty Components of an Incomplete Ensemble of Climate Projections Using Data Augmentation. Journal of Climate. <doi:10.1175/JCLI-D-18-0606.1>.

Examples

##########################################################################
# SYNTHETIC SCENARIOS
##########################################################################
# create nS=3 fictive climate scenarios with 2 GCMs and 2 RCMs, for a period of nY=20 years
n=20
t=1:n/n

# GCM effects (sums to 0 for each t)
effGCM1 = t*2
effGCM2 = t*-2

# RCM effects (sums to 0 for each t)
effRCM1 = t*1
effRCM2 = t*-1

# These climate scenarios are a sum of effects and a random gaussian noise
scenGCM1RCM1 = effGCM1 + effRCM1 + rnorm(n=n,sd=0.5)
scenGCM1RCM2 = effGCM1 + effRCM2 + rnorm(n=n,sd=0.5)
scenGCM2RCM1 = effGCM2 + effRCM1 + rnorm(n=n,sd=0.5)
Y.synth = rbind(scenGCM1RCM1,scenGCM1RCM2,scenGCM2RCM1)

# Here, scenAvail indicates that the first scenario is obtained with the combination of the
# GCM "GCM1" and RCM "RCM1", the second scenario is obtained with the combination of
# the GCM "GCM1" and RCM "RCM2" and the third scenario is obtained with the combination
# of the GCM "GCM2" and RCM "RCM1".
scenAvail.synth = data.frame(GCM=c('GCM1','GCM1','GCM2'),RCM=c('RCM1','RCM2','RCM1'))

##########################################################################
# RUN QUALYPSO
##########################################################################
# call main QUALYPSO function: two arguments are mandatory:
# - Y: Climate projections for nS scenarios and nY time steps. if Y is a matrix nS x nY, we
# run QUALYPSO nY times, for each time step. If Y is an array nG x nS x nY, for nG grid points,
# we run QUALYPSO nG times, for each grid point, for one time step specified using the argument
# iFut
# - scenAvail: matrix or data.frame of available combinations nS x nEff. The number of
# characteristics nEff corresponds to the number of main effects that will be included in the
# ANOVA model. In the following example, we have nEff=2 main effects corresponding to the GCMs
# and RCMs.

# Many options can be specified in the argument "listOption". When ANOVAmethod=="QUALYPSO"
# a Bayesian inference is performed. Here, we change the default values for nBurn and nKeep
# in order to speed up computation time for this small example. However, it must be noticed
# that convergence and sampling of the posterior distributions often require higher values
#  for these two arguments.
listOption = list(nBurn=100,nKeep=100,ANOVAmethod="QUALYPSO",quantilePosterior=c(0.025,0.5,0.975))

# run QUALYPSO
QUALYPSO.synth = QUALYPSO(Y=Y.synth, scenAvail=scenAvail.synth, X=2001:2020, listOption=listOption)

##########################################################################
# SOME PLOTS
##########################################################################
# plot grand mean
plotQUALYPSOgrandmean(QUALYPSO.synth,xlab="Years")

# plot main GCM effects
plotQUALYPSOeffect(QUALYPSO.synth,nameEff="GCM",xlab="Years")

# plot main RCM effects
plotQUALYPSOeffect(QUALYPSO.synth,nameEff="RCM",xlab="Years")

# plot fraction of total variance for the differences sources of uncertainty
plotQUALYPSOTotalVarianceDecomposition(QUALYPSO.synth,xlab="Years")

# plot mean prediction and total variance with the differences sources of uncertainty
plotQUALYPSOMeanChangeAndUncertainties(QUALYPSO.synth,xlab="Years")

#____________________________________________________________
# EXAMPLE OF QUALYPSO WHEN THE PREDICTOR IS TIME
#____________________________________________________________

# list of options
listOption = list(typeChangeVariable='abs')

# call QUALYPSO
QUALYPSO.time = QUALYPSO(Y=Y,scenAvail=scenAvail,X=X_time_vec,
                         Xfut=Xfut_time,listOption=listOption)

# grand mean effect
plotQUALYPSOgrandmean(QUALYPSO.time,xlab="Years")

# main GCM effects
plotQUALYPSOeffect(QUALYPSO.time,nameEff="GCM",xlab="Years")

# main RCM effects
plotQUALYPSOeffect(QUALYPSO.time,nameEff="RCM",xlab="Years")

# mean change and associated uncertainties
plotQUALYPSOMeanChangeAndUncertainties(QUALYPSO.time,xlab="Years")

# variance decomposition
plotQUALYPSOTotalVarianceDecomposition(QUALYPSO.time,xlab="Years")

#____________________________________________________________
# EXAMPLE OF QUALYPSO WHEN THE PREDICTOR IS THE GLOBAL TEMPERATURE
#____________________________________________________________

# list of options
listOption = list(typeChangeVariable='abs')

# call QUALYPSO
QUALYPSO.globaltas = QUALYPSO(Y=Y,scenAvail=scenAvail,X=X_globaltas,
                              Xfut=Xfut_globaltas,listOption=listOption)

# grand mean effect
plotQUALYPSOgrandmean(QUALYPSO.globaltas,xlab="Global warming (Celsius)")

# main GCM effects
plotQUALYPSOeffect(QUALYPSO.globaltas,nameEff="GCM",xlab="Global warming (Celsius)")

# main RCM effects
plotQUALYPSOeffect(QUALYPSO.globaltas,nameEff="RCM",xlab="Global warming (Celsius)")

# mean change and associated uncertainties
plotQUALYPSOMeanChangeAndUncertainties(QUALYPSO.globaltas,xlab="Global warming (Celsius)")

# variance decomposition
plotQUALYPSOTotalVarianceDecomposition(QUALYPSO.globaltas,xlab="Global warming (Celsius)")


QUALYPSO documentation built on Oct. 24, 2023, 9:07 a.m.