diff --git a/R/multiscaleSVDxpts.R b/R/multiscaleSVDxpts.R index 23457961..ac092a71 100644 --- a/R/multiscaleSVDxpts.R +++ b/R/multiscaleSVDxpts.R @@ -1,6 +1,8 @@ #' No fail SVD function that switches to rsvd if svd fails #' #' This function performs SVD on a matrix using the built-in svd function in R. +#' The matrix will be divided by its maximum value before computing the SVD for +#' the purposes of numerical stability (optional). #' If svd fails, it automatically switches to random svd from the rsvd package. #' svd may fail to converge when the matrix condition number is high; this can #' be checked with the kappa function. @@ -8,6 +10,7 @@ #' @param x Matrix to perform SVD on #' @param nu Number of left singular vectors to return (default: min(nrow(x), ncol(x))) #' @param nv Number of right singular vectors to return (default: min(nrow(x), ncol(x))) +#' @param dividebymax boolean #' #' @return A list containing the SVD decomposition of x #' @@ -16,14 +19,14 @@ #' nc <- 10 #' u <- ba_svd( avgU, nu = nc, nv = 0)$u #' @export -ba_svd <- function(x, nu = min(nrow(x), ncol(x)), nv = min(nrow(x), ncol(x))) { +ba_svd <- function(x, nu = min(nrow(x), ncol(x)), nv = min(nrow(x), ncol(x)), dividebymax=FALSE ) { tryCatch( expr = { - svd(x, nu = nu, nv = nv) + svd(x/max(x), nu = nu, nv = nv) }, error = function(e) { message("svd failed, using rsvd instead") - rsvd(x, nu = nu, nv = nv) + rsvd(x/max(x), nu = nu, nv = nv) } ) } @@ -4271,8 +4274,8 @@ simlr.perm <- function(voxmats, smoothingMatrices, iterations = 10, sparsenessQu #' @export rvcoef <- function(X, Y) { # Normalize matrices - X_norm <- scale(X, center = TRUE, scale = TRUE) - Y_norm <- scale(Y, center = TRUE, scale = TRUE) + X_norm <- scale(X/max(X), center = TRUE, scale = TRUE) + Y_norm <- scale(Y/max(Y), center = TRUE, scale = TRUE) # Compute cross-product matrix C <- t(X_norm) %*% Y_norm diff --git a/man/ba_svd.Rd b/man/ba_svd.Rd index 66dab3e5..f6d0c706 100644 --- a/man/ba_svd.Rd +++ b/man/ba_svd.Rd @@ -4,7 +4,12 @@ \alias{ba_svd} \title{No fail SVD function that switches to rsvd if svd fails} \usage{ -ba_svd(x, nu = min(nrow(x), ncol(x)), nv = min(nrow(x), ncol(x))) +ba_svd( + x, + nu = min(nrow(x), ncol(x)), + nv = min(nrow(x), ncol(x)), + dividebymax = FALSE +) } \arguments{ \item{x}{Matrix to perform SVD on} @@ -12,12 +17,16 @@ ba_svd(x, nu = min(nrow(x), ncol(x)), nv = min(nrow(x), ncol(x))) \item{nu}{Number of left singular vectors to return (default: min(nrow(x), ncol(x)))} \item{nv}{Number of right singular vectors to return (default: min(nrow(x), ncol(x)))} + +\item{dividebymax}{boolean} } \value{ A list containing the SVD decomposition of x } \description{ This function performs SVD on a matrix using the built-in svd function in R. +The matrix will be divided by its maximum value before computing the SVD for +the purposes of numerical stability (optional). If svd fails, it automatically switches to random svd from the rsvd package. svd may fail to converge when the matrix condition number is high; this can be checked with the kappa function.