gevdfit.R 25.9 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
# This file contains the functions:
# - gev.d.fit, gev.d.init for fitting
# - gev.d.diag for diagnostic plots
# - gev.d.params, gev.d.rl for calculation of parameters/ return levels
# and the documentation of the example data

#### gev.d.fit ####

#' @title Maximum-likelihood Fitting of the duration dependent GEV Distribution
#' @description Modified \code{\link[ismev]{gev.fit}} function for Maximum-likelihood fitting 
#' for the duration dependent generalized extreme 
#' value distribution, following Koutsoyiannis et al. (1988), including generalized linear 
#' modelling of each parameter.
#' @param xdat A vector containing maxima for different durations. 
#' This can be obtained from \code{\link{IDF.agg}}.
#' @param ds A vector of aggregation levels corresponding to the maxima in xdat.
#' @param ydat A matrix of covariates for generalized linear modelling of the parameters 
#' (or NULL (the default) for stationary fitting). The number of rows should be the same as the 
#' length of xdat.
#' @param  mul,sigl,shl,thetal,etal Numeric vectors of integers, giving the columns of ydat that contain
#'  covariates for generalized linear modelling of the parameters (or NULL (the default) 
#'  if the corresponding parameter is stationary).
#'  Parameters are: modified location, scale_0, shape, duration offset, duration exponent repectively.
#' @param mulink,siglink,shlink,thetalink,etalink Inverse link functions for generalized linear 
#' modelling of the parameters.
#' @param muinit,siginit,shinit,thetainit,etainit initial values as numeric of length 
#' equal to total number of parameters 
#' used to model the parameters. If NULL (the default) is given, initial parameters are obtained 
#' internally by fitting the GEV seperately for each duration and applying a linear model to optain the 
#' duration dependency of the location and shape parameter.
#' @param show Logical; if TRUE (the default), print details of the fit.
#' @param method The optimization method used in \code{\link{optim}}.
#' @param maxit The maximum number of iterations.
#' @param ... Other control parameters for the optimization.
#' @return A list containing the following components. 
#' A subset of these components are printed after the fit. 
#' If show is TRUE, then assuming that successful convergence is indicated, 
#' the components nllh, mle and se are always printed. 
#' \item{nllh}{single numeric giving the negative log-likelihood value.} 
#' \item{mle}{numeric vector giving the MLE's for the modified location, scale_0, shape, 
#' duration offset and duration exponent, resp.} 
#' \item{se}{numeric vector giving the standard errors for the MLE's (in the same order).}
#' \item{trans}{An logical indicator for a non-stationary fit.}
#' \item{model}{A list with components mul, sigl, shl, thetal and etal.}
#' \item{link}{A character vector giving inverse link functions.}
#' \item{conv}{The convergence code, taken from the list returned by \code{\link{optim}}. 
#' A zero indicates successful convergence.}
#' \item{data}{data is standardized to standart Gumbel.} 
#' \item{cov}{The covariance matrix.} 
#' @seealso \code{\link{dgev.d}}, \code{\link{IDF.agg}}, \code{\link{gev.fit}}, \code{\link{optim}}
#' @export
#' @importFrom stats optim 
#' 
#' @examples 
#' # sampled random data from d-gev with covariates
#' # GEV parameters:
#' # mu = 4 + 0.2*cov1 +0.5*cov2
#' # sigma = 2+0.5*cov1
#' # xi = 0.5
#' # theta = 0
#' # eta = 0.5
#' 
#' data('example',package ='IDF')
#' 
#' gev.d.fit(xdat=example$dat,ds = example$d,ydat=as.matrix(example[,c('cov1','cov2')])
#' ,mul=c(1,2),sigl=1)

gev.d.fit<-
  function(xdat, ds, ydat = NULL, mul = NULL, sigl = NULL, shl = NULL, thetal = NULL, etal = NULL, 
           mulink = identity, siglink = identity, shlink = identity, thetalink = identity, etalink = identity,  
           muinit = NULL, siginit = NULL, shinit = NULL, thetainit = NULL, etainit = NULL,
           show = TRUE, method = "Nelder-Mead", maxit = 10000, ...)
  {
    #
    # obtains mles etc for d-gev distn
    #
    
    # test for NA values:
    if(any(is.na(xdat))) stop('xdat contains NA values. NA values need to be removed first.')
    
    z <- list()
    # number of parameters (betas) to estimate for each parameter: 
    npmu <- length(mul) + 1
    npsc <- length(sigl) + 1
    npsh <- length(shl) + 1
    npth <- length(thetal) + 1
    npet <- length(etal) + 1
    z$trans <- FALSE  # indicates if fit is non-stationary
    
    # calculate initial values for mu.d, sigma_0, xi, eta using IDF.init:  (thetainit=0)
    init.vals <- gev.d.init(xdat,ds,ifelse(is.null(thetainit),0,thetainit[1]))
Jana Ulrich's avatar
Jana Ulrich committed
92
93
94
    # TODO: transform initial values with link function
    # use make.link()
    
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
    
    # generate covariates matrices: 
    if (is.null(mul)) {
      mumat <- as.matrix(rep(1, length(xdat)))
      if (is.null(muinit)) 
        muinit <- init.vals$mu
    }else {
      z$trans <- TRUE
      mumat <- cbind(rep(1, length(xdat)), ydat[, mul])
      if (is.null(muinit)) 
        muinit <- c(init.vals$mu, rep(0, length(mul)))
    }
    if (is.null(sigl)) {
      sigmat <- as.matrix(rep(1, length(xdat)))
      if (is.null(siginit)) 
        siginit <- init.vals$sigma
    }else {
      z$trans <- TRUE
      sigmat <- cbind(rep(1, length(xdat)), ydat[, sigl])
      if (is.null(siginit)) 
        siginit <- c(init.vals$sigma, rep(0, length(sigl)))
    }
    if (is.null(shl)) {
      shmat <- as.matrix(rep(1, length(xdat)))
      if (is.null(shinit)) 
        shinit <- init.vals$xi 
    }else {
      z$trans <- TRUE
      shmat <- cbind(rep(1, length(xdat)), ydat[, shl])
      if (is.null(shinit)) 
        shinit <- c(init.vals$xi, rep(0, length(shl)))
    }
    if (is.null(thetal)) {
      thmat <- as.matrix(rep(1, length(xdat)))
      if (is.null(thetainit))  
        thetainit <- 0
    }else {
      z$trans <- TRUE
      thmat <- cbind(rep(1, length(xdat)), ydat[, thetal])
      if (is.null(thetainit))  
        thetainit <- c(0, rep(0, length(thetal)))
    }
    if (is.null(etal)) {
      etmat <- as.matrix(rep(1, length(xdat)))
      if (is.null(etainit)) 
        etainit <- init.vals$eta
    }else {
      z$trans <- TRUE
      etmat <- cbind(rep(1, length(xdat)), ydat[, etal])
      if (is.null(etainit)) 
        etainit <- c(init.vals$eta, rep(0, length(etal)))
    }
    
    z$model <- list(mul, sigl, shl, thetal, etal)
    z$link <- deparse(substitute(c(mulink, siglink, shlink, thetalink, etalink)))
Jana Ulrich's avatar
Jana Ulrich committed
150
    # TODO: store as functions not as strings!
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
    init <- c(muinit, siginit, shinit, thetainit, etainit)
    
    # function to calculate neg log-likelihood:
    gev.lik <- function(a) {
      # computes neg log lik of d-gev model
      mu <- mulink(mumat %*% (a[1:npmu]))
      sigma <- siglink(sigmat %*% (a[seq(npmu + 1, length = npsc)]))
      xi <- shlink(shmat %*% (a[seq(npmu + npsc + 1, length = npsh)]))
      theta <- thetalink(thmat %*% (a[seq(npmu + npsc + npsh + 1, length = npth)]))
      eta <- etalink(etmat %*% (a[seq(npmu + npsc + npsh + npth + 1, length = npet)]))
      
      ds.t <- ds+theta
      sigma.d <- sigma/(ds.t^eta)
      y <- xdat/sigma.d - mu
      y <- 1 + xi * y
      
      if(any(eta <= 0) ||any(theta <= -0.01) || any(sigma.d <= 0) || any(y <= 0)) return(10^6)
      sum(log(sigma.d)) + sum(y^(-1/xi)) + sum(log(y) * (1/xi + 1))
    }
    
    #####################################################################################
    # derivations of nll after d-gev-parameters (for boosting):
    
174
    # get parameters from covariates (and a vector containing predictors)
175
176
177
178
179
180
    # mu <- mulink(mumat %*% (a[1:npmu]))
    # sigma <- siglink(sigmat %*% (a[seq(npmu + 1, length = npsc)]))
    # xi <- shlink(shmat %*% (a[seq(npmu + npsc + 1, length = npsh)]))
    # theta <- thetalink(thmat %*% (a[seq(npmu + npsc + npsh + 1, length = npth)]))
    # eta <- etalink(etmat %*% (a[seq(npmu + npsc + npsh + npth + 1, length = npet)]))
    # 
Jana Ulrich's avatar
Jana Ulrich committed
181
182
183
184
185
186
187
188
189
    # xsd <- xdat*(ds+theta)^eta/sigma
    # y <- 1 + xi * (xsd - mu)
    # 
    # nll <- log(sigma.d) + y^(-1/xi) + log(y) * (1/xi + 1)
    # dnll.mu <- y^(-1/xi-1)-(1+xi)/y
    # dnll.sigma <- xsd*y^(-1/xi-1)/sigma -(1+xi)*xsd/sigma/y+1/sigma
    # dnll.xi <- y^(-1/xi)*(log(y)/xi^2-(xsd-mu)/(xi*y))-log(y)/xi^2+(1/xi+1)*(xsd-mu)/y
    # dnll.theta <-  (-eta*xsd*y^(-1/xi-1)+eta*(1+xi)*xsd/y-eta)/(ds+theta)
    # dnll.eta <- (-xsd*y^(-1/xi-1)+(1+xi)*xsd/y-1)*log(ds+theta)
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
    #####################################################################################
    
    
    # finding minimum of log-likelihood:
    x <- optim(init, gev.lik, hessian = TRUE, method = method,
               control = list(maxit = maxit, ...))
    
    # saving output parameters:
    z$conv <- x$convergence
    mut <- mulink(mumat %*% (x$par[1:npmu]))
    sc0 <- siglink(sigmat %*% (x$par[seq(npmu + 1, length = npsc)]))
    xi <- shlink(shmat %*% (x$par[seq(npmu + npsc + 1, length = npsh)]))
    theta <- thetalink(thmat %*% (x$par[seq(npmu + npsc + npsh + 1, length = npth)]))
    eta <- etalink(etmat %*% (x$par[seq(npmu + npsc + npsh + npth + 1, length = npet)]))
    z$nllh <- x$value
    # normalize data to standart gumbel:
    sc.d <- sc0/((ds+theta)^eta)
    z$data <-  - log(as.vector((1 + xi * (xdat/sc.d-mut))^(-1/xi))) 
    z$mle <- x$par
209
    test <- try(              # catch error 
210
    z$cov <- solve(x$hessian) # invert hessian to get estimation on var-covar-matrix
211
212
213
214
215
    ,silent = TRUE )
    if("try-error" %in% class(test)){
      warning("Hessian could not be inverted. NAs were produced.")
      z$cov <- matrix(NA,length(z$mle),length(z$mle))
        }
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
    z$se <- sqrt(diag(z$cov)) # sqrt(digonal entries) = standart error of mle's 
    z$vals <- cbind(mut, sc0, xi, theta, eta)
    colnames(z$vals) <- c('mut','sigma0','xi','theta','eta')
    z$ds <- ds
    if(show) {
      if(z$trans) # for nonstationary fit
        print(z[c(2, 3, 4)]) # print model, link, conv
      else print(z[4]) # for stationary fit print only conv
      if(!z$conv) # if fit converged 
        print(z[c(5, 7, 9)]) # print nll, mle, se
    }
    class( z) <- "gev.d.fit"
    invisible(z)
}


#### gev.d.init ####

# function to get initial values for gev.d.fit:
# obtain initial values 
# by fitting every duration seperately

# possible ways to improve:
# take given initial values into accout, if there are any
# xi -> mean vs. median ... how do we improve that?
# mu_tilde -> is not very good for small sample sizes yet
# improved inital value for eta, by fitting both mu~d and sigma~d in log-log scale

#' @title get initial values for gev.d.fit
#' @description obtain initial values by fitting every duration seperately
#' @param xdat vector of maxima for differnt durations
#' @param ds vector of durations belonging to maxima in xdat
#' @param thetainit initial parameter for theta
#' @return list of initail values for mu_tilde, sigma_0, xi, eta
#' @importFrom stats lm 
#' @importFrom ismev gev.fit
#' @keywords internal 

gev.d.init <- function(xdat,ds,thetainit){
  
  durs <- unique(ds)
  mles <- matrix(NA, nrow=length(durs), ncol= 3)
  for(i in 1:length(durs)){
Jana Ulrich's avatar
Jana Ulrich committed
259
260
    test <- try(mles[i,] <- gev.fit(xdat[ds==durs[i]],show = FALSE)$mle,silent = TRUE)
    if("try-error" %in% class(test)){mles[i,] <- rep(NA,3)}
261
262
263
264
265
266
267
268
269
270
271
272
  }
  # get values for sig0 and eta (also mu_0) from linear model in log-log scale
  lmsig <- lm(log(mles[,2])~log(durs+thetainit))
  lmmu <- lm(log(mles[,1])~log(durs+thetainit))
  
  # sig0 <- exp Intercept
  siginit <- exp(lmsig$coefficients[[1]])
  # eta <- mean of negativ slopes 
  etainit <- mean(c(-lmsig$coefficients[[2]],-lmmu$coefficients[[2]]))
  # mean of mu_d/sig_d 
  # could try:
  # mu0/sig0 is also an estimate but needs to be weighted in mean
Jana Ulrich's avatar
Jana Ulrich committed
273
  muinit <- mean(c(mles[,1]/mles[,2]),na.rm = TRUE) #exp(lmmu$coefficients[[1]])/exp(lmsig$coefficients[[1]])
274
  # mean of shape parameters 
Jana Ulrich's avatar
Jana Ulrich committed
275
  shinit <- mean(mles[,3],na.rm = TRUE)
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
  
  return(list(mu=muinit,sigma=siginit,xi=shinit,eta=etainit))
}


#### gev.d.diag ####

#' Diagnostic Plots for d-gev Models
#'
#' @description  Produces diagnostic plots for d-gev models using 
#' the output of the function \code{\link{gev.d.fit}}. Values for different durations can be plotted in 
#' different colors of with different symbols.
#' @param fit object returned by \code{\link{gev.d.fit}}
#' @param subset an optional vector specifying a subset of observations to be used in the plot
#' @param cols optional either one value or vector of same length as \code{unique(durations)} to
#' specify the colors of plotting points. 
#' The default uses the \code{rainbow} function.
#' @param pch optional either one value or vector of same length as \code{unique(durations)} containing
#' integers or symbols to specify the plotting points.
#' @param which string containing 'both', 'pp' or 'qq' to specify, which plots should be produced.
#' @param mfrow vector specifying layout of plots. If both plots should be produced seperately,
#' set to \code{c(1,1)}.
#' @param legend logical indicating if legends should be plotted
Jana Ulrich's avatar
Jana Ulrich committed
299
300
301
#' @param title character vector of length 2, giving the titles for the pp- and the qq-plot
#' @param emp.lab,mod.lab character string containing names for empirical and model axis
#' @param ... additional parameters passed on to the plotting function 
302
303
304
305
306
307
308
309
310
311
312
313
314
315
#'
#' @export
#' @importFrom graphics plot abline par title
#' @importFrom grDevices rainbow
#'
#' @examples
#' data('example',package ='IDF')
#' 
#' fit <- gev.d.fit(xdat=example$dat,ds = example$d,ydat=as.matrix(example[,c('cov1','cov2')])
#'                  ,mul=c(1,2),sigl=1)
#' # diagnostic plots for complete data                
#' gev.d.diag(fit)    
#' # diagnostic plots for subset of data (e.g. one station)            
#' gev.d.diag(fit,subset = example$cov1==1)
Jana Ulrich's avatar
Jana Ulrich committed
316
317
318
gev.d.diag <- function(fit,subset=NULL,cols=NULL,pch=NULL,which='both',mfrow=c(1,2),legend=TRUE,
                       title=c('Residual Probability Plot','Residual Quantile Plot'),
                       emp.lab='Empirical',mod.lab='Model',...){
319
  # check parameter:
Jana Ulrich's avatar
Jana Ulrich committed
320
  if(!is.element(which,c('both','pp','qq'))) stop("Parameter 'which'= ",which,
321
322
323
324
                                                 " but only 'both','pp' or 'qq' are allowed.")
  # subset data
  df <- data.frame(data=fit$data,ds=fit$ds)
  if(!is.null(subset))df <- df[subset,]
Jana Ulrich's avatar
Jana Ulrich committed
325
326
327
328
329
  # get single durations
  durs <- sort(unique(df$ds))
  # rescale durations to assign colors
  df$cval <- sapply(df$ds,function(d){which(durs==d)})

330
331
332
333
334
335
336
337
338
339
  # sort data 
  df <- df[order(df$data),]
  
  # plotting position
  n <- length(df$data)
  px <- (1:n)/(n + 1)

  # create plots:
  if(which=='both') par(mfrow=mfrow) # 2 subplots
  # colors and symbols
Jana Ulrich's avatar
Jana Ulrich committed
340
  if(is.null(cols))cols <- rainbow(length(durs))
341
342
343
344
345
  if(is.null(pch))pch <- df$cval
  
  if(which=='both'|which=='pp'){
    # pp
    plot(px, exp( - exp( - df$data)), xlab =
Jana Ulrich's avatar
Jana Ulrich committed
346
347
348
           emp.lab, ylab = mod.lab,col=cols[df$cval],pch=pch,...)
    abline(0, 1, col = 1,lwd=1)
    title(title[1])
Jana Ulrich's avatar
Jana Ulrich committed
349
    if(legend){legend('bottomright',legend = round(durs,digits = 2),pch=pch,
Jana Ulrich's avatar
Jana Ulrich committed
350
                      col = cols[1:length(durs)],title = 'Durations',ncol = 2)}
351
352
353
354
  }
  if(which=='both'|which=='qq'){
    # qq
    plot( - log( - log(px)), df$data, ylab =
Jana Ulrich's avatar
Jana Ulrich committed
355
356
357
            emp.lab, xlab = mod.lab,col=cols[df$cval],pch=pch,...)
    abline(0, 1, col = 1,lwd=1)
    title(title[2])
Jana Ulrich's avatar
Jana Ulrich committed
358
    if(legend){legend('bottomright',legend = round(durs,digits = 2),pch=pch,
Jana Ulrich's avatar
Jana Ulrich committed
359
                      col = cols[1:length(durs)],title = 'Durations',ncol = 2)}
360
361
362
363
364
365
366
367
368
369
370
  }
  if(which=='both') par(mfrow=c(1,1)) # reset par
}

#### gev.d.params ####

#' Calculate gev(d) parameters from \code{gev.d.fit} output
#'
#' @description function to calculate mut, sigma0, xi, theta, eta 
#' (modified location, scale, shape, duration offset, duration exponent) 
#' from results of \code{\link{gev.d.fit}} with covariates
371
#' @param fit fit object returned by \code{gev.d.fit} or \code{gev.fit}
372
#' @param ydat A matrix containing the covariates in the same order as used in \code{gev.d.fit}.
373
#' @seealso \code{\link{dgev.d}}
374
#' @return data.frame containing mu_tilde, sigma0, xi, theta, eta (or mu, sigma, xi for gev.fit objects)
375
376
377
378
379
380
#' @export
#' 
#' @examples
#' data('example',package = 'IDF')
#' fit <- gev.d.fit(example$dat,example$d,ydat = as.matrix(example[,c("cov1","cov2")])
#'                  ,mul = c(1,2),sigl = 1)
381
382
383
384
#' gev.d.params(fit = fit,ydat = cbind(c(0.9,1),c(0.5,1)))


gev.d.params <- function(fit,ydat){
385
386
  if(!class(fit)%in%c("gev.d.fit","gev.fit"))stop("'fit' must be an object returned by 'gev.d.fit' or 'gev.fit'.")
  if(!is.matrix(ydat))stop("'ydat' must be of class matrix.")
Jana Ulrich's avatar
Jana Ulrich committed
387
388
389
  if(!fit$trans){warning('No vglm for parameters. Max. likelihood estimates are returned.')
    return(data.frame(mut=rep(fit$mle[1],dim(ydat)[1]),sig0=fit$mle[2],xi=fit$mle[3]
                      ,theta=fit$mle[4],eta=fit$mle[5]))}
390
391
392
  n.par <- max(sapply(fit$model,function(x){return(ifelse(is.null(x),0,max(x)))}))
  if(n.par>ncol(ydat))stop("Covariates-Matrix 'ydat' has ",ncol(ydat), " columns, but ", n.par," are required.")
  
393
  #ydat <- rbind(0,ydat) # no error in case ncols=1
394
  
395
396
397
398
399
400
  # number of parameters
  npmu <- length(fit$model[[1]]) + 1
  npsc <- length(fit$model[[2]]) + 1
  npsh <- length(fit$model[[3]]) + 1
  if(class(fit)=="gev.d.fit"){npth <- length(fit$model[[4]]) + 1}
  if(class(fit)=="gev.d.fit"){npet <- length(fit$model[[5]]) + 1}
401
  
402
403
404
405
406
407
  # link functions
  mulink <- eval(parse(text=fit$link))[[1]]
  siglink <- eval(parse(text=fit$link))[[2]]
  shlink <- eval(parse(text=fit$link))[[3]]
  if(class(fit)=="gev.d.fit"){thetalink <- eval(parse(text=fit$link))[[4]]}
  if(class(fit)=="gev.d.fit"){etalink <- eval(parse(text=fit$link))[[5]]}
408
  
409
410
411
412
413
414
  # covariates matrices
  mumat <- cbind(rep(1, dim(ydat)[1]), matrix(ydat[, fit$model[[1]]],dim(ydat)[1],npmu-1))
  sigmat <- cbind(rep(1, dim(ydat)[1]), matrix(ydat[, fit$model[[2]]],dim(ydat)[1],npsc-1))
  shmat <- cbind(rep(1, dim(ydat)[1]), matrix(ydat[, fit$model[[3]]],dim(ydat)[1],npsh-1))
  if(class(fit)=="gev.d.fit"){thmat <- cbind(rep(1, dim(ydat)[1]), matrix(ydat[, fit$model[[4]]],dim(ydat)[1],npth-1))}
  if(class(fit)=="gev.d.fit"){etmat <- cbind(rep(1, dim(ydat)[1]), matrix(ydat[, fit$model[[5]]],dim(ydat)[1],npet-1))}
415
  
416
417
418
419
420
421
  # calculate parameters
  mut <- mulink(mumat %*% (fit$mle[1:npmu]))
  sc0 <- siglink(sigmat %*% (fit$mle[seq(npmu + 1, length = npsc)]))
  xi <- shlink(shmat %*% (fit$mle[seq(npmu + npsc + 1, length = npsh)]))
  if(class(fit)=="gev.d.fit"){theta <- thetalink(thmat %*% (fit$mle[seq(npmu + npsc + npsh + 1, length = npth)]))}
  if(class(fit)=="gev.d.fit"){eta <- etalink(etmat %*% (fit$mle[seq(npmu + npsc + npsh + npth + 1, length = npet)]))}
422
  
423
424
  if(class(fit)=="gev.d.fit"){return(data.frame(mut=mut,sig0=sc0,xi=xi,theta=theta,eta=eta))
  }else{return(data.frame(mu=mut,sig=sc0,xi=xi))}
425
426
427
}


428

429
430
#### gev.d.rl ####

Jana Ulrich's avatar
Jana Ulrich committed
431
#' Calculate Returnlevel 
432
#'
Jana Ulrich's avatar
Jana Ulrich committed
433
#' @description calculate Returnlevel for chosen duration and return period 
434
435
436
437
438
439
440
#' from \code{\link{gev.d.fit}} parameters
#' @param params list of parameters mu_tilde, sigma0, xi, theta, eta (single values and/or compatible matrices)
#' as obtained from \code{\link{gev.d.fit}} or \code{\link{gev.d.params}}
#' @param p.d numeric vector of length 2 containing one value the for exeedance probability p = 1-1/RP 
#' and one value for the duration at which to calculate the return level
#'
#' @return one return level value or matrix with return levels (depending on input to params)
441
#' unit: e.g. mm/h 
442
443
444
445
446
447
448
#' @export
#' 
#' @examples 
#' data('example',package = 'IDF')
#' fit <- gev.d.fit(example$dat,example$d,ydat = as.matrix(example[,c("cov1","cov2")])
#'                  ,mul = c(1,2),sigl = 1)
#' ### calculate rl on grid:
449
450
451
#' # covariates values
#' cov1 <- rep(seq(-1,1,0.1),11)
#' cov2 <- rep(seq(0,1,0.1),each=21)
452
#' # calculate parameters of d-gev on grid, based on output of gev.d.fit
453
#' par <- gev.d.params(fit = fit,ydat = cbind(cov1,cov2))
454
455
#' # calculate 100 year (p=0.99) rl for a duration of d=24 hours
#' rl <- gev.d.rl(params = par,p.d = c(0.99,24))
456
457
#' dim(rl) <- c(21,11)
#' # rl map:
458
459
460
461
462
463
#' image(x=seq(-1,1,0.1),y=seq(0,1,0.1),z=rl,xlab = 'cov1', ylab = 'cov2')
gev.d.rl <- function(params,p.d){
  sigma.d <- params[[2]]/((p.d[2]+params[[4]])^params[[5]])
  mu <- params[[1]]*sigma.d
  yt <- -1/log(p.d[1])
  rl <- mu+sigma.d/params[[3]]*(yt^params[[3]]-1)
464
  return(rl)
465
466
}

Jana Ulrich's avatar
Jana Ulrich committed
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
### alternativ function that can also calculate delta-method-sd:

#' Calculate Returnlevel with standart deviation (through delta method)
#'
#'@description calculate returnlevel for chosen duration and return period 
#' for \code{\link{gev.d.fit}} objects
#' @param p single numeric value containing the non-exceedance probability p=1-1/T
#' @param fit \code{\link{gev.d.fit}} object
#' @param ds numeric vector of durations for which to calculate the return levels
#' @param ydat matrix containing the covariates in the same order as used in \code{\link{gev.d.fit}}
#' @param sd logical, if TRUE, standart deviation of the returnlevels obtained with the delta method
#' are also returned
#'
#' @return either vector or matrix containing the p-quantile or list containing both the p-quantile and its standart deviation
#' for durations in ds (and different rows for the covariates in ydat) 
#' @export
#'
#' @examples
#' data('example',package = 'IDF')
#' fit <- gev.d.fit(example$dat,example$d,ydat = as.matrix(example[,c("cov1","cov2")])
#'                  ,mul = c(1,2),sigl = 1)
#' # calculate 100 year (p=0.99) rl for a duration of d=24 hours 
#' rl <- gev.d.rl2(p = 0.99,fit = fit,ds = 1:30,ydat = matrix(0.7,1,2),sd=TRUE)
#' plot(1:30,rl$q,log='xy',type = 'l',xlab = 'Duration',ylab = 'Intensity',ylim = c(5,100))
#' # add 0.95-confidence intervals
#' lines(1:30,rl$q-1.96*rl$sd,lty=2)
#' lines(1:30,rl$q+1.96*rl$sd,lty=2)
gev.d.rl2<-function(p,fit,ds,ydat,sd=FALSE){
  if(length(p)>1){
    warning("'p' can only be a single value. Only first element of 'p' is used.")
    p <- p[1]
  }
  # parameter estimates
  if(!fit$trans){     # model without covariates
    mut <- fit$mle[1]
    sig0 <- fit$mle[2]
    xi <- fit$mle[3]
    theta <- fit$mle[4]
    eta <- fit$mle[5]
  }else{                    # model with covariates
    n.cov <- max(sapply(fit$model,function(x){return(ifelse(is.null(x),0,max(x)))}))
    if(is.null(ydat)){      # if no ydat, 0's will be used as values for covariates 
      warning("No covariates-matrix was given. Zeros are used as values of covariates.")
      ydat <- matrix(0,1,n.cov)
    }
    if(!is.matrix(ydat)){
      warning("'ydat' was converted to a matrix.")
      ydat <- matrix(ydat,1,n.cov)
    }
    yrows <- nrow(ydat)
    ds <- matrix(ds,nrow(ydat),length(ds),byrow = TRUE)
    params <- gev.d.params(fit,ydat)  # calculate parameters
    
    ydat <- cbind(rep(1,yrows),ydat) # add first column of ones
    
    mut <- params$mut
    sig0 <- params$sig0
    xi <- params$xi
    theta <- params$theta
    eta <- params$eta
  }
  # quantile
  y.xi <- (-1/log(p))^xi
  dtheta <- ds+theta
  sigma <- sig0*dtheta^(-eta)
  q <- mut*sigma+sigma/xi*(y.xi-1)
  if(fit$trans){colnames(q) <- ds[1,]}
  if(!sd){return(q)} # only qantiles
  # partial derivatives
  dq.dmut <- sigma
  dq.dsig0 <- mut*dtheta^(-eta)-dtheta^(-eta)/xi*(y.xi-1)
  dq.dxi <- -sigma/xi^2*(y.xi-1)-sigma/xi*y.xi*log(-log(p))
  dq.dtheta <- mut*sig0*dtheta^(-eta-1)-sig0*(-eta*dtheta^(-eta-1))/xi*(y.xi-1)
  dq.deta <- -mut*sigma*log(dtheta)+sigma/xi*log(dtheta)*(y.xi-1)
  if(!fit$trans){
    dq.dpar <- rbind(dq.dmut,dq.dsig0,dq.dxi,dq.dtheta,dq.deta) 
  }else{
    dq.dpar <- lapply(1:yrows,function(i){
      rbind(dq.dmut[i,],dq.dsig0[i,],dq.dxi[i,],dq.dtheta[i,],dq.deta[i,])})
  }
  # calculate sd
  if(!fit$trans){     # model without covariates 
    sd<-apply(dq.dpar,2,function(d.vec) sqrt(t(d.vec) %*% fit$cov %*% d.vec))
  }else{                    # model with covariates
    n.par <- sapply(fit$model,length)+1
    cov.ind <- sapply(fit$model,function(vec){c(1,vec+1)})
    cov.ind <- do.call(c,cov.ind)               # index of ydat
    
    calc.sd <- function(j){
      inner <- ydat[j,cov.ind]                  # inner derivative
      outer <- dq.dpar[[j]]
      outer <- sapply(1:5,function(i){          # outer derivative
        matrix(rep(outer[i,],n.par[i]),n.par[i],length(ds[1,]),byrow = TRUE)})
      outer <- do.call(rbind,outer)
      sd<-apply(outer,2,function(d.vec) sqrt(t(d.vec*inner) %*% fit$cov %*% (d.vec*inner)))
    }
    sd <- lapply(1:yrows,calc.sd)               # calculate sd for every row in ydat
    sd <- do.call(rbind,sd)
    colnames(sd) <- ds[1,]
  }
  return(list(q=q,sd=sd))     # qantiles and standart deviation
}
569

570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605

#### gev.d2stdgumbel ####

#' Transform data to standart gumbel
#'
#' @param xdat A vector containing maxima for different durations.
#' @param ds A vector of aggregation levels corresponding to the maxima in xdat.
#' @param params list of parameters mu_tilde, sigma0, xi, theta, eta 
#' as obtained from \code{\link{gev.d.params}}
#'
#' @return Vector containing transformed data.
#' @export
#'
#' @examples 
#' data('example',package = 'IDF')
#' # fit subset
#' ind <- sample(1:10, length(example$dat), replace=TRUE)
#' fit.subs <- gev.d.fit(example$dat[ind!=1],example$d[ind!=1]
#'                       ,ydat = as.matrix(example[ind!=1,c("cov1","cov2")])
#'                       ,mul = c(1,2),sigl = 1)
#' # calculate parameters for unfitted values
#' par <- gev.d.params(fit = fit.subs
#'                       ,ydat = as.matrix(example[ind==1,c("cov1","cov2")]))
#' # transform unfitted values to standart gumbel
#' sg.data <- gev.d2stdgumbel(xdat = example$dat[ind==1]
#'                       ,ds = example$d[ind==1],params = par)
#' # check unfitted values agains standart gumbel
#' gev.d.diag(data.frame(data=sg.data,ds=example$d[ind==1]),pch=20)
gev.d2stdgumbel <- function(xdat,ds,params){
  sc.d <- params$sig0/((ds+params$theta)^params$eta)
  sg.data <-  - log(as.vector((1 + params$xi * (xdat/sc.d-params$mut))^(-1/params$xi))) 
  return(sg.data)
}



606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
#### example data ####

#' Sampled data for duration dependent GEV
#'
#' A dataset containing:
#' \itemize{
#'   \item \code{$xdat}: 'annual' maxima values
#'   \item \code{$ds}: corresponding durations
#'   \item \code{$cov1}, \code{$cov2}: covariates}
#' GEV parameters:
#' \itemize{
#'   \item mu = 4 + 0.2*cov1 +0.5*cov2
#'   \item sigma = 2+0.5*cov1
#'   \item xi = 0.5
#'   \item theta = 0
#'   \item eta = 0.5}
#'
#' @docType data
#' @keywords datasets
#' @name example
#' @usage data('example',package ='IDF')
NULL