covls.R 10.2 KB
Newer Older
1
#' @include procmod.R
2
#' @include procmod_frame.R
3
#' @include multivariate.R
4 5
#' @import doParallel
#' @import foreach
6
#' @import stats
7 8 9 10 11
#'
#' @author Christelle Gonindard-Melodelima
#' @author Eric Coissac
NULL

12 13 14
library(doParallel)
registerDoParallel(1)

15 16
#' Compute the trace of a square matrix.
#'
17 18 19
#' The trace of a square matrix is defined as the sum
#' of its diagonal elements.
#'
20 21 22 23
#' @param X a square matrix
#' @return the trace of X
#'
#' @examples
24 25
#' m <- matrix(1:16, nrow = 4)
#' ProcMod:::.Trace(m)
26 27
#' @note Internal function do not use.
#'
28
#' @rdname internal.Trace
29 30 31
#' @author Eric Coissac
#' @author Christelle Gonindard-Melodelima
#'
32
.Trace <- function(X) sum(diag(X))
33

34 35 36 37 38 39 40 41 42 43 44
#' Convert a covariance matrix to a correlation matric
#'
#' @param c a square covariance/variance matrix
#' @return correlation matrix corresponding to \code{c}
#'
#' @note Internal function do not use.
#'
#' @rdname internal.var2cor
#' @author Eric Coissac
#' @author Christelle Gonindard-Melodelima
#'
45 46 47 48 49 50
.var2cor <- function(c) {
  v <- sqrt(diag(c))
  vv <- v %o% v
  c / vv
}

Eric Coissac committed
51 52 53 54 55
#' Procrustean Correlation, and Variance / Covariance Matrices.
#'
#' \code{varls}, \code{corls}, \code{corls.partial} compute the procrustean
#' variance / covariance, correlation, or partial correlation matrices
#' between a set of real matrices and \code{\link[stats]{dist}} objects.
56
#'
57
#' Procrustean covariance between two matrices X and Y, is defined as the sum
Eric Coissac committed
58 59 60 61 62 63 64 65
#' of the singular values of the X'Y matrix \insertCite{Gower:71:00,Lingoes:74:00}{ProcMod}.
#' Both the X and Y matrices must have the same number of rows.
#'
#' The variances and covariances and correlations are corrected
#' to avoid over fitting \insertCite{Coissac-Eric:19:00}{ProcMod}.
#'
#' Partial correlation coefficients are computed by inverting the correlation followed
#' by a normalisation by the diagonal of the inverted matrix.
66
#'
Eric Coissac committed
67 68 69
#' The inputs must be numeric matrices or \code{\link[stats]{dist}} object.
#' The set of input matrices can be aggregated un a
#' \code{\link[ProcMod]{procmod.frame}}.
70
#'
Eric Coissac committed
71
#' Before computing the coefficients, matrices are projected into an
72 73
#' orthogonal space using the \code{\link[ProcMod]{ortho}} function.
#'
Eric Coissac committed
74 75 76 77 78 79 80 81 82
#' The denominator n - 1 is used which gives an unbiased estimator of the
#' (co)variance for i.i.d. observations.
#'
#' Scaling a covariance matrix into a correlation one can be achieved in many ways,
#' mathematically most appealing by multiplication with a diagonal matrix from left
#' and right, or more efficiently by using sweep(.., FUN = "/") twice.
#' The \code{\link[stats]{cov2cor}} function is even a bit more efficient,
#' and provided mostly for didactical reasons.
#'
83 84 85 86 87 88 89 90 91 92
#' @references{
#'  \insertRef{Gower:71:00}{ProcMod}
#'
#'  \insertRef{Lingoes:74:00}{ProcMod}
#'
#'  \insertRef{Coissac-Eric:19:00}{ProcMod}
#' }
#'
#' @param ...   the set of matrices or a \code{\link[ProcMod]{procmod.frame}}
#'              object.
93 94
#' @param nrand number of randomisation used to estimate the mean
#'              covariance observed between two random matrix.
95 96 97
#'              If rand is \code{NULL} or equal to \code{0}, no correction
#'              is estimated and the raw procrustean covariances are
#'              estimated.
98
#' @param p.adjust.method the multiple test correction method used
99
#'              to adjust p values. \code{p.adjust.method}
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
#'              belongsone of the folowing values: \code{"holm"},
#'              \code{"hochberg"}, \code{"hommel"}, \code{"bonferroni"},
#'              \code{"BH"}, \code{"BY"}, \code{"fdr"},  \code{"none"}.
#'              The default is,set to \code{"holm"}.
#'
#' @return a \code{procmod.varls} object which corresponds to a numeric
#'              matrix annotated by several attributes.
#'
#'              The following attribute is always added:
#'
#'              - \code{nrand} an integer value indicating the number of
#'                randomisations used to estimate the mean of the random
#'                covariance.
#'
#'              When \code{nrand} is greater than 0 a couple of attributes
#'              is added:
#'
#'              - \code{rcovls} a numeric matrix containing the estimation
#'                of the mean of the random covariance.
#'
#'              - \code{p.value} a numeric matrix containing the estimations
#'                of the p.values of tests checking that the observed
#'                covariance is larger than the mean of the random covariance.
#'                p.values are corrected for multiple tests according to the
#'                method specified by the \code{p.adjust.method} parameter.
125 126
#'
#' @examples
127 128 129 130 131
#' # Build Three matrices of 3 rows.
#' A <- matrix(1:9, nrow = 3)
#' B <- matrix(10:15, nrow = 3)
#' C <- matrix(20:31, nrow = 3)
#' # compute the variance covariance matrix
132 133
#' varls(A, B, C)
#' varls(A = A, B = B, C = C)
134 135
#' data = procmod.frame(A = A, B = B, C = C)
#' varls(data)
Eric Coissac committed
136
#'
137 138
#' @author Eric Coissac
#' @author Christelle Gonindard-Melodelima
Eric Coissac committed
139
#'
140 141
#' @seealso \code{\link[stats]{p.adjust}}
#'
Eric Coissac committed
142 143 144 145 146
#' @rdname varls
#' @name varls
#' @aliases varls
#' @aliases corls
#' @aliases corls.partial
147
#' @export
148 149 150 151 152
varls <- function(...,
                  nrand = 100,
                  p.adjust.method = "holm") {
  if (!is.null(nrand) && length(nrand) > 1) {
    stop("nrand must be a single numeric value or the NULL value")
153
  }
154

155 156 157
  if (nrand == 0) nrand <- NULL

  xs <- list(...)
158

159 160
  if (length(xs) == 1) {
    x <- xs[[1]]
161
    if (is.procmod.frame(x)) {
162
      xs <- x
163
    } else {
164
      xs <- procmod.frame(x)
165 166 167
    }
  }
  else {
168
    xs <- as.procmod.frame(xs)
169 170
  }

171
  x_names <- names(xs)
172

173
  xs <- ortho(xs)
174

175 176 177
  nx <- length(xs)
  n <- nrow(xs)
  n1 <- n - 1
178

179
  cov_xxs <- Matrix::Matrix(0, nrow = nx, ncol = nx)
180 181

  # Computes the standard Covls covariance matrix
182
  for (i in seq_len(nx))
183 184
    for (j in i:nx)
      cov_xxs[i, j] <- sum(svd(crossprod(xs[[i]],xs[[j]]))$d) / n1
185

186
  cov_xxs <- as.matrix(Matrix::forceSymmetric(cov_xxs, uplo = "U"))
187

188 189 190
  # Computes Covls under null hypothesis
  if (!is.null(nrand)) {
    n_r_greater <- array(0, dim = dim(cov_xxs))
191

192 193 194 195
    v_xs <- vector(mode = "list", nx)
    for (i in seq_len(nx))
      v_xs[[i]] <- var(xs[[i]])

Eric Coissac committed
196 197
    if (! getDoParRegistered()) registerDoParallel(1)

198 199 200
    s_cov_xxs <- foreach(k = seq_len(nrand),
                         .combine = cbind) %dopar% {
      s1_cov_xxs <- matrix(0, nrow = nx, ncol = nx)
201
      r_xs <- vector(mode = "list", nx)
202 203 204 205 206 207 208 209 210 211 212
      r_ys <- vector(mode = "list", nx)
      for (i in seq_len(nx)) {
        r_xs[[i]] <- MASS::mvrnorm(n,
                                   rep(0,nrow(v_xs[[i]])),
                                   Sigma = v_xs[[i]],
                                   empirical = TRUE)
        r_ys[[i]] <- MASS::mvrnorm(n,
                                   rep(0,nrow(v_xs[[i]])),
                                   Sigma = v_xs[[i]],
                                   empirical = TRUE)
      }
213 214 215

      for (i in seq_len(nx))
        for (j in i:nx)
216
          s1_cov_xxs[i, j] <- sum(svd(crossprod(r_xs[[i]],r_ys[[j]]))$d)
217

218
      s1_cov_xxs  / n1
219

220
                         }
221

222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238
    dim(s_cov_xxs) = c(nx, nx, nrand)

    s_cov_xxs <- s_cov_xxs

    r_cov_xxs <- apply(s_cov_xxs,
      MARGIN = c(1,2),
      FUN = mean
    )

    for (k in seq_len(nrand))
      n_r_greater <- n_r_greater + (s_cov_xxs[, , k] >= cov_xxs)

    r_cov_xxs <- as.matrix(Matrix::forceSymmetric(r_cov_xxs, uplo = "U"))
    cov_xxs <- pmax(cov_xxs - r_cov_xxs, 0)

    colnames(r_cov_xxs) <- x_names
    rownames(r_cov_xxs) <- x_names
239

240
    p_values <- n_r_greater / nrand
241

242 243 244 245 246 247 248 249 250 251 252 253 254
    c_p_values <- p.adjust(p_values[upper.tri(p_values, diag = FALSE)],
      method = p.adjust.method
    )

    p_values[upper.tri(p_values, diag = FALSE)] <- c_p_values
    p_values <- as.matrix(Matrix::forceSymmetric(p_values, uplo = "U"))

    colnames(p_values) <- x_names
    rownames(p_values) <- x_names

    attr(cov_xxs, "rcovls") <- as.matrix(r_cov_xxs)
    attr(cov_xxs, "p.value") <- p_values
    attr(cov_xxs, "nrand") <- nrand
255
  }
256

257 258
  colnames(cov_xxs) <- x_names
  rownames(cov_xxs) <- x_names
259

260
  make_subS3Class(cov_xxs, "procmod.varls")
261 262 263
}


Eric Coissac committed
264
#' @rdname varls
265
#' @export
266 267 268 269 270 271
corls <- function(..., nrand = 100,
                  p.adjust.method = "holm") {
  cov <- varls(...,
    nrand = nrand,
    p.adjust.method = p.adjust.method
  )
272
  s <- sqrt(diag(cov))
273 274
  vv <- s %o% s
  rls <- cov / vv
Eric Coissac committed
275
  class(rls) <- "matrix"
276 277


278 279
  if (!is.null(attr(cov, "rcovls"))) {
    attr(rls, "nrand") <- attr(cov, "nrand")
280
  }
281

282
  make_subS3Class(rls, "procmod.corls")
283 284
}

Eric Coissac committed
285
#' @rdname varls
286
#' @export
287 288
corls.partial <- function(..., nrand = 100) {
  rls <- corls(..., nrand = nrand)
289
  C <- solve(rls)
290 291 292
  S <- sqrt(diag(C))
  D <- S %o% S
  rp <- C / D
293

294 295
  attr(rp, "nrand") <- attr(rls, "nrand")
  attr(rp, "rcovls") <- attr(rls, "rcovls")
296

297
  make_subS3Class(rp, "procmod.corls")
298 299
}

Eric Coissac committed
300
#' Print procrustean variance / covariance matrix.
301
#'
302 303 304 305 306 307
#' @param x a \code{procmod.varls}
#'          object
#' @param ... other parameters passed to other functions
#'
#' @seealso \code{\link[ProcMod]{varls}}
#'
308 309 310
#' @author Eric Coissac
#' @author Christelle Gonindard-Melodelima
#' @export
311 312
print.procmod.varls <- function(x, ...) {
  class(x) <- "matrix"
313
  attr(x, "nrand") <- NULL
314
  attr(x, "rcovls") <- NULL
315
  attr(x, "p.value") <- NULL
316 317 318
  print(x)
}

319 320 321 322 323 324 325 326 327 328 329
#' @author Eric Coissac
#' @author Christelle Gonindard-Melodelima
#' @export
`$.procmod.varls` <- function(x, name) {
  attr(x,name)
}

#' @author Eric Coissac
#' @author Christelle Gonindard-Melodelima
#' @export
names.procmod.varls <- function(x) {
Eric Coissac committed
330 331 332 333
  n <- names(attributes(x))
  bn <- grep(pattern = "^(dim|class)",
             x = n)
  n[-bn]
334 335 336
}


Eric Coissac committed
337
#' Print procrustean correlation matrix.
338
#'
339 340 341 342 343 344
#' @param x a \code{procmod.corls}
#'          object
#' @param ... other parameters passed to other functions
#'
#' @seealso \code{\link[ProcMod]{corls}}
#'
345 346 347
#' @author Eric Coissac
#' @author Christelle Gonindard-Melodelima
#' @export
348 349
print.procmod.corls <- function(x, ...) {
  class(x) <- "matrix"
350
  attr(x, "nrand") <- NULL
351
  attr(x, "rcovls") <- NULL
352 353
  attr(x, "rcorls") <- NULL
  attr(x, "p.value") <- NULL
354 355
  print(x)
}
356 357 358 359 360 361 362 363

#' @author Eric Coissac
#' @author Christelle Gonindard-Melodelima
#' @export
`$.procmod.corls` <- function(x, name) {
  attr(x,name)
}

Eric Coissac committed
364 365 366 367 368 369 370 371 372 373 374 375

#' @author Eric Coissac
#' @author Christelle Gonindard-Melodelima
#' @export
names.procmod.corls <- function(x) {
  n <- names(attributes(x))
  bn <- grep(pattern = "^(dim|class)",
             x = n)
  n[-bn]
}