2021-11-12 17:22:45 +00:00
|
|
|
solve.gep <- function(A, B, d = nrow(A)) {
|
|
|
|
isrB <- matpow(B, -0.5)
|
|
|
|
|
|
|
|
if (requireNamespace("RSpectra", quietly = TRUE) && d < nrow(A)) {
|
|
|
|
eig <- RSpectra::eigs_sym(isrB %*% A %*% isrB, d)
|
|
|
|
} else {
|
|
|
|
eig <- eigen(isrB %*% A %*% isrB, symmetric = TRUE)
|
|
|
|
}
|
|
|
|
|
|
|
|
list(vectors = isrB %*% eig$vectors, values = eig$values)
|
|
|
|
}
|
|
|
|
|
|
|
|
POI.lambda.max <- function(A, d = 1L, method = c('POI-C', 'POI-L', 'FastPOI-C', 'FastPOI-L')) {
|
|
|
|
method <- match.arg(method)
|
|
|
|
|
|
|
|
if (method %in% c('POI-C', 'POI-L')) {
|
|
|
|
A2 <- apply(apply(A^2, 2, sort, decreasing = TRUE), 2, cumsum)
|
|
|
|
lambda.max <- sqrt(apply(A2[1:d, , drop = FALSE], 1, max))
|
|
|
|
if (method == 'POI-C') {
|
|
|
|
lambda.max[d]
|
|
|
|
} else {
|
|
|
|
lambda.max[1]
|
|
|
|
}
|
|
|
|
} else {
|
|
|
|
if (requireNamespace("RSpectra", quietly = TRUE)) {
|
|
|
|
vec <- RSpectra::eigs_sym(A, d)$vectors
|
|
|
|
} else {
|
|
|
|
vec <- eigen(A, symmetric = TRUE)$vectors
|
|
|
|
}
|
|
|
|
if (method == 'FastPOI-C') {
|
|
|
|
sqrt(max(rowSums(vec^2)))
|
|
|
|
} else { # 'FastPOI-L'
|
|
|
|
max(abs(vec))
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
2020-06-10 14:35:27 +00:00
|
|
|
#' Penalysed Orthogonal Iteration.
|
|
|
|
#'
|
|
|
|
#' @param lambda Default: 0.75 * lambda_max for FastPOI-C method.
|
|
|
|
#'
|
|
|
|
#' @note use.C required 'poi.so' beeing dynamicaly loaded.
|
|
|
|
#' dyn.load('../tensor_predictors/poi.so')
|
2021-10-29 16:16:40 +00:00
|
|
|
#'
|
|
|
|
#' @export
|
2021-11-12 17:22:45 +00:00
|
|
|
POI <- function(A, B, d = 1L, sparsity = 0.5,
|
|
|
|
method = c('POI-C', 'FastPOI-C'), # TODO: Maybe implement the the lasso loss too
|
|
|
|
iter.outer = 100L, iter.inner = 500L,
|
|
|
|
tol = sqrt(.Machine$double.eps)
|
|
|
|
) {
|
|
|
|
method <- match.arg(method)
|
|
|
|
|
|
|
|
# Compute penalty parameter lambda
|
|
|
|
lambda <- sparsity * POI.lambda.max(A, d, method)
|
|
|
|
|
|
|
|
# Ensure RHS B to be positive definite
|
|
|
|
if (missing(B)) {
|
|
|
|
B <- diag(nrow(A))
|
|
|
|
} else {
|
|
|
|
rankB <- qr(B, tol)$rank
|
|
|
|
if (rankB < nrow(B)) {
|
|
|
|
diag(B) <- diag(B) + log(nrow(B)) / rankB
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
# In case of zero penalty compute ordinary GEP solution
|
|
|
|
if (lambda == 0) {
|
|
|
|
eig <- solve.eig(A, B, d)
|
|
|
|
return(structure(list(
|
|
|
|
U = eig$vectors, d = eig$values,
|
|
|
|
lambda = 0, call = match.call()
|
|
|
|
), class = c("tensor_predictor", "POI")))
|
|
|
|
}
|
|
|
|
|
|
|
|
# Set initial values
|
|
|
|
if (requireNamespace("RSpectra", quietly = TRUE) && nrow(A) < d) {
|
2021-11-04 12:05:15 +00:00
|
|
|
Delta <- RSpectra::eigs_sym(A, d)$vectors
|
2020-06-10 14:35:27 +00:00
|
|
|
} else {
|
2021-11-04 12:05:15 +00:00
|
|
|
Delta <- eigen(A, symmetric = TRUE)$vectors[, 1:d, drop = FALSE]
|
2020-06-10 14:35:27 +00:00
|
|
|
}
|
2021-11-12 17:22:45 +00:00
|
|
|
Q <- Delta
|
|
|
|
Z <- matrix(0, nrow(Q), ncol(Q))
|
2020-06-10 14:35:27 +00:00
|
|
|
|
2021-11-12 17:22:45 +00:00
|
|
|
# Outer loop (iteration)
|
|
|
|
for (i in seq_len(iter.outer)) {
|
|
|
|
Q.last <- Q # for break condition
|
2020-06-10 14:35:27 +00:00
|
|
|
|
2021-11-12 17:22:45 +00:00
|
|
|
# Step 1: Solve B Z_i = A Q_{i-1} for Z_i
|
|
|
|
Delta <- crossprod(A, Q)
|
|
|
|
# Inner Loop
|
|
|
|
for (j in seq_len(iter.inner)) {
|
|
|
|
Z.last <- Z # for break condition
|
|
|
|
|
|
|
|
traces <- Delta - B %*% Z + diag(B) * Z
|
|
|
|
Z <- traces * (pmax(1 - lambda / sqrt(rowSums(traces^2)), 0) / diag(B))
|
|
|
|
|
|
|
|
# Inner break condition
|
|
|
|
if (norm(Z.last - Z, 'F') < tol) {
|
|
|
|
break
|
2020-06-10 14:35:27 +00:00
|
|
|
}
|
|
|
|
}
|
|
|
|
|
2021-11-12 17:22:45 +00:00
|
|
|
# Step 2: QR decomposition of Z_i = Q_i R_i.
|
|
|
|
if (d == 1L) {
|
|
|
|
Z.norm <- norm(Z, 'F')
|
|
|
|
if (Z.norm < tol) {
|
|
|
|
Q <- matrix(0, p, d)
|
|
|
|
} else {
|
|
|
|
Q <- Z / Z.norm
|
|
|
|
}
|
2020-06-10 14:35:27 +00:00
|
|
|
} else {
|
2021-11-12 17:22:45 +00:00
|
|
|
# Detect zero columns.
|
|
|
|
zero.col <- colSums(abs(Z)) < tol
|
|
|
|
if (all(zero.col)) {
|
|
|
|
Q <- matrix(0, p, d)
|
|
|
|
} else if (any(zero.col)) {
|
|
|
|
Q <- matrix(0, p, d)
|
|
|
|
Q[, !zero.col] <- qr.Q(qr(Z[, !zero.col]))
|
|
|
|
} else {
|
|
|
|
Q <- qr.Q(qr(Z))
|
|
|
|
}
|
2020-06-10 14:35:27 +00:00
|
|
|
}
|
2021-11-12 17:22:45 +00:00
|
|
|
|
|
|
|
# In case of fast POI, only one iteration
|
|
|
|
if (startsWith(method, 'Fast')) {
|
|
|
|
break
|
|
|
|
}
|
|
|
|
if (norm(tcrossprod(Q, Q) - tcrossprod(Q.last, Q.last), 'F') < tol) {
|
|
|
|
break
|
2020-06-10 14:35:27 +00:00
|
|
|
}
|
|
|
|
}
|
|
|
|
|
2021-11-12 17:22:45 +00:00
|
|
|
# TODO: Finish with transformation to original solution U of
|
|
|
|
# A U = B U Lambda.
|
|
|
|
|
|
|
|
structure(list(
|
|
|
|
Z = Z, Q = Q,
|
|
|
|
lambda = lambda,
|
|
|
|
call = match.call()
|
|
|
|
), class = c("tensor_predictor", "POI"))
|
2020-06-10 14:35:27 +00:00
|
|
|
}
|
2021-11-12 17:22:45 +00:00
|
|
|
|
|
|
|
|
|
|
|
# POI.bak <- function(A, B, d,
|
|
|
|
# lambda = 0.75 * sqrt(max(rowSums(Delta^2))),
|
|
|
|
# update.tol = 1e-3,
|
|
|
|
# tol = 100 * .Machine$double.eps,
|
|
|
|
# maxit = 400L,
|
|
|
|
# # maxit.outer = maxit,
|
|
|
|
# maxit.inner = maxit,
|
|
|
|
# use.C = FALSE,
|
|
|
|
# method = 'FastPOI-C') {
|
|
|
|
|
|
|
|
# # TODO:
|
|
|
|
# stopifnot(method == 'FastPOI-C')
|
|
|
|
|
|
|
|
# if (requireNamespace("RSpectra", quietly = TRUE)) {
|
|
|
|
# Delta <- RSpectra::eigs_sym(A, d)$vectors
|
|
|
|
# } else {
|
|
|
|
# Delta <- eigen(A, symmetric = TRUE)$vectors[, 1:d, drop = FALSE]
|
|
|
|
# }
|
|
|
|
|
|
|
|
# # Set initial value.
|
|
|
|
# Z <- Delta
|
|
|
|
|
|
|
|
# # Step 1: Optimization.
|
|
|
|
# # The "inner" optimization loop, aka repeated coordinate optimization.
|
|
|
|
# if (use.C) {
|
|
|
|
# Z <- .Call('FastPOI_C_sub', A, B, Delta, lambda, as.integer(maxit.inner),
|
|
|
|
# PACKAGE = 'tensorPredictors')
|
|
|
|
# } else {
|
|
|
|
# p <- nrow(Z)
|
|
|
|
# for (iter.inner in 1:maxit.inner) {
|
|
|
|
# Zold <- Z
|
|
|
|
# for (g in 1:p) {
|
|
|
|
# a <- Delta[g, ] - B[g, ] %*% Z + B[g, g] * Z[g, ]
|
|
|
|
# a_norm <- sqrt(sum(a^2))
|
|
|
|
# if (a_norm > lambda) {
|
|
|
|
# Z[g, ] <- a * ((1 - lambda / a_norm) / B[g, g])
|
|
|
|
# } else {
|
|
|
|
# Z[g, ] <- 0
|
|
|
|
# }
|
|
|
|
# }
|
|
|
|
# if (norm(Z - Zold, 'F') < update.tol) {
|
|
|
|
# break
|
|
|
|
# }
|
|
|
|
# }
|
|
|
|
# }
|
|
|
|
|
|
|
|
# # Step 2: QR decomposition.
|
|
|
|
# if (d == 1L) {
|
|
|
|
# Z_norm <- sqrt(sum(Z^2))
|
|
|
|
# if (Z_norm < tol) {
|
|
|
|
# Q <- matrix(0, p, d)
|
|
|
|
# } else {
|
|
|
|
# Q <- Z / Z_norm
|
|
|
|
# }
|
|
|
|
# } else {
|
|
|
|
# # Detect zero columns.
|
|
|
|
# zeroColumn <- colSums(abs(Z)) < tol
|
|
|
|
# if (all(zeroColumn)) {
|
|
|
|
# Q <- matrix(0, p, d)
|
|
|
|
# } else if (any(zeroColumn)) {
|
|
|
|
# Q <- matrix(0, p, d)
|
|
|
|
# Q[, !zeroColumn] <- qr.Q(qr(Z))
|
|
|
|
# } else {
|
|
|
|
# Q <- qr.Q(qr(Z))
|
|
|
|
# }
|
|
|
|
# }
|
|
|
|
|
|
|
|
# list(Z = Z, Q = Q, iter.inner = if (use.C) NA else iter.inner,
|
|
|
|
# lambda = lambda)
|
|
|
|
# }
|