Skip to content
Open
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
3 changes: 2 additions & 1 deletion DESCRIPTION
Original file line number Diff line number Diff line change
Expand Up @@ -5,7 +5,8 @@ Date: 2021-07-13
Authors@R: c(
person("Patrick", "Breheny", role=c("aut","cre"), email="patrick-breheny@uiowa.edu", comment=c(ORCID="0000-0002-0650-1119")),
person("Yaohui", "Zeng", role="ctb"),
person("Ryan", "Kurth", role="ctb"))
person("Ryan", "Kurth", role="ctb"),
person("Jamie", "Wallis", role="ctb"))
Depends: R (>= 3.1.0)
Imports: Matrix
Suggests: knitr, rmarkdown, splines, survival, tinytest
Expand Down
2 changes: 1 addition & 1 deletion R/auc.R
Original file line number Diff line number Diff line change
Expand Up @@ -4,7 +4,7 @@ AUC.cv.grpsurv <- function(obj, ...) {
stop("The 'survival' package is needed for AUC() to work. Please install it.", call. = FALSE)
}
if (packageVersion("survival") < "3.2.10") stop("AUC.cv.ncvsurv requires version 3.2.10 of 'survival' package or higher", call.=FALSE)
S <- survival::Surv(obj$fit$time, obj$fit$fail)
S <- survival::Surv(obj$fit$start_time, obj$fit$stop_time, obj$fit$fail)
res <- apply(obj$Y, 2, concord, y = S)
num <- res['concordant',] + 0.5*res['tied.x',] + 0.5*res['tied.y',] + 0.5*res['tied.xy',]
num/sum(res[,1])
Expand Down
2 changes: 1 addition & 1 deletion R/cv-grpsurv.R
Original file line number Diff line number Diff line change
Expand Up @@ -11,7 +11,7 @@ cv.grpsurv <- function(X, y, group, ..., nfolds=10, seed, fold, se=c('quick', 'b

# Get standardized X, y
X <- fit$XG$X
y <- cbind(fit$time, fit$fail)
y <- cbind(fit$start_time, fit$stop_time, fit$fail)
returnX <- list(...)$returnX
if (is.null(returnX) || !returnX) fit$X <- NULL

Expand Down
12 changes: 7 additions & 5 deletions R/grpsurv.R
Original file line number Diff line number Diff line change
Expand Up @@ -42,24 +42,25 @@ grpsurv <- function(X, y, group=1:ncol(X), penalty=c("grLasso", "grMCP", "grSCAD

# Set up lambda
if (missing(lambda)) {
lambda <- setupLambdaCox(XG$X, Y$time, Y$fail, XG$g, penalty, alpha, lambda.min, nlambda, XG$m)
lambda <- setupLambdaCox(XG$X, Y, XG$g, penalty, alpha, lambda.min, nlambda, XG$m)
user.lambda <- FALSE
} else {
nlambda <- length(lambda)
user.lambda <- TRUE
}

## Fit
n <- length(Y$time)
n <- length(Y$fail)
p <- ncol(XG$X)
K <- as.integer(table(XG$g))
K0 <- as.integer(if (min(XG$g)==0) K[1] else 0)
K1 <- as.integer(if (min(XG$g)==0) cumsum(K) else c(0, cumsum(K)))
start_order = as.numeric(order(Y$start_time))
if (bilevel) {
res <- .Call("lcdfit_cox", XG$X, Y$fail, penalty, K1, K0, lambda, alpha, eps, 0, gamma, tau, as.integer(max.iter),
res <- .Call("lcdfit_cox", XG$X, Y$start_time, start_order, Y$stop_time, Y$fail, penalty, K1, K0, lambda, alpha, eps, 0, gamma, tau, as.integer(max.iter),
XG$m, as.integer(dfmax), as.integer(gmax), as.integer(warn), as.integer(user.lambda))
} else {
res <- .Call("gdfit_cox", XG$X, Y$fail, penalty, K1, K0, lambda, alpha, eps, as.integer(max.iter),
res <- .Call("gdfit_cox", XG$X, Y$start_time, start_order, Y$stop_time, Y$fail, penalty, K1, K0, lambda, alpha, eps, as.integer(max.iter),
as.double(gamma), XG$m, as.integer(dfmax), as.integer(gmax), as.integer(warn), as.integer(user.lambda))
}
b <- matrix(res[[1]], p, nlambda)
Expand Down Expand Up @@ -101,7 +102,8 @@ grpsurv <- function(X, y, group=1:ncol(X), penalty=c("grLasso", "grMCP", "grSCAD
df = df,
iter = iter,
group.multiplier = XG$m,
time = Y$time,
start_time = Y$start_time,
stop_time = Y$stop_time,
fail = Y$fail,
order = Y$ind,
linear.predictors = sweep(Eta, 2, colMeans(Eta), '-')),
Expand Down
18 changes: 16 additions & 2 deletions R/loss.R
Original file line number Diff line number Diff line change
Expand Up @@ -22,15 +22,29 @@ loss.grpreg <- function(y, yhat, family) {
val
}
loss.grpsurv <- function(y, eta, total=TRUE) {
ind <- order(y[,1])
d <- as.integer(y[ind,2])
ind <- order(y[,2])
d <- as.integer(y[ind,3])
if (is.matrix(eta)) {
eta <- eta[ind, , drop=FALSE]
r <- apply(eta, 2, function(x) rev(cumsum(rev(exp(x)))))
} else {
eta <- as.matrix(eta[ind])
r <- as.matrix(rev(cumsum(rev(exp(eta)))))
}
n = nrow(y)
current_sum = 0
start_idx = n
stop_idx = n
start_o = order(y[,1])
while(start_idx>0 && stop_idx >0){
if(y[start_o[start_idx],1]<y[ind[stop_idx],2]){
r[stop_idx,] = r[stop_idx,] - current_sum
stop_idx = stop_idx - 1
}else{
current_sum = current_sum + exp(eta[start_o[start_idx],])
start_idx = start_idx - 1
}
}
if (total) {
return(-2*(crossprod(d, eta) - crossprod(d, log(r))))
} else {
Expand Down
19 changes: 15 additions & 4 deletions R/newS.R
Original file line number Diff line number Diff line change
@@ -1,6 +1,17 @@
newS <- function(y) {
ind <- order(y[,1])
list(time = as.double(y[ind,1]),
fail = as.double(y[ind,2]),
ind = ind)
if(ncol(y) == 2){
ind <- order(y[,1])
list(start_time = as.double(rep(0,nrow(y))),
stop_time = as.double(y[ind,1]),
fail = as.double(y[ind,2]),
ind = ind)
}else if(ncol(y)==3){
ind <- order(y[,2])
list(start_time = as.double(y[ind,1]),
stop_time = as.double(y[ind,2]),
fail = as.double(y[ind,3]),
ind = ind)
}else{
stop('y must be specified as Surv(time,fail) or Surv(start,stop,fail)')
}
}
4 changes: 2 additions & 2 deletions R/predict-surv.R
Original file line number Diff line number Diff line change
Expand Up @@ -97,7 +97,7 @@ predict.grpsurv <- function(object, X,
a <- ifelse(object$fail, (1-w/r)^(1/w), 1)
S0 <- c(1, cumprod(a))
H0 <- c(0, cumsum(1-a))
x <- c(0, object$time)
x <- c(0, object$stop_time)
for (i in 1:nrow(eta)) {
S <- S0^exp(eta[i,j])
H <- H0*exp(eta[i,j])
Expand All @@ -113,7 +113,7 @@ predict.grpsurv <- function(object, X,
if (type %in% c('survival', 'hazard')) {
if (nrow(eta)==1) val <- val[[1]]
class(val) <- c('grpsurv.func', class(val))
attr(val, 'time') <- object$time
attr(val, 'time') <- object$stop_time
} else if (type == 'median') {
val <- drop(val)
}
Expand Down
2 changes: 1 addition & 1 deletion R/residuals.R
Original file line number Diff line number Diff line change
Expand Up @@ -23,7 +23,7 @@ residuals.grpreg <- function(object, lambda, which=1:length(object$lambda), drop
# Calculate matrix of residuals
if (inherits(object, 'grpsurv')) {
for (j in 1:length(object$lambda)) {
h <- suppressWarnings(predict(object, which=j, type='hazard')(object$time))
h <- suppressWarnings(predict(object, which=j, type='hazard')(object$stop_time))
M <- object$fail - h * exp(object$linear.predictors)
R <- sign(M) * sqrt(-2*(M + object$fail*log(object$fail-M)))
R[h==0,] <- 0
Expand Down
21 changes: 17 additions & 4 deletions R/setupLambdaCox.R
Original file line number Diff line number Diff line change
@@ -1,4 +1,4 @@
setupLambdaCox <- function(X, y, Delta, group, penalty, alpha, lambda.min, nlambda, group.multiplier) {
setupLambdaCox <- function(X, y, group, penalty, alpha, lambda.min, nlambda, group.multiplier) {
n <- nrow(X)
p <- ncol(X)

Expand All @@ -8,13 +8,26 @@ setupLambdaCox <- function(X, y, Delta, group, penalty, alpha, lambda.min, nlamb
if (K1[1]!=0) {
SURV <- get("Surv", asNamespace("survival"))
COXPH <- get("coxph", asNamespace("survival"))
nullFit <- COXPH(SURV(y, Delta) ~ X[, group==0, drop=FALSE])
nullFit <- COXPH(SURV(y$start_time, y$stop_time, y$fail) ~ X[, group==0, drop=FALSE])
eta <- nullFit$linear.predictors
rsk <- rev(cumsum(rev(exp(eta))))
s <- Delta - exp(eta)*cumsum(Delta/rsk)
current_sum = 0
start_idx = n
stop_idx = n
start_o = order(y[,1])
while(start_idx>0 && stop_idx >0){
if(y[start_o[start_idx],1]<y[ind[stop_idx],2]){
rsk[stop_idx] = rsk[stop_idx] - current_sum
stop_idx = stop_idx - 1
}else{
current_sum = current_sum + exp(eta[start_o[start_idx]])
start_idx = start_idx - 1
}
}
s <- y$fail - exp(eta)*cumsum(y$fail/rsk)
} else {
w <- 1/(n-(1:n)+1)
s <- Delta - cumsum(Delta*w)
s <- y$fail - cumsum(y$fail*w)
}

## Determine lambda.max
Expand Down
20 changes: 19 additions & 1 deletion src/gdfit_cox.c
Original file line number Diff line number Diff line change
Expand Up @@ -45,7 +45,7 @@ void gd_cox(double *b, double *x, double *r, double *eta, double v, int g,
Free(z);
}

SEXP gdfit_cox(SEXP X_, SEXP d_, SEXP penalty_, SEXP K1_, SEXP K0_, SEXP lambda, SEXP alpha_, SEXP eps_, SEXP max_iter_, SEXP gamma_, SEXP group_multiplier, SEXP dfmax_, SEXP gmax_, SEXP warn_, SEXP user_) {
SEXP gdfit_cox(SEXP X_, SEXP start_, SEXP start_o_, SEXP stop_, SEXP d_, SEXP penalty_, SEXP K1_, SEXP K0_, SEXP lambda, SEXP alpha_, SEXP eps_, SEXP max_iter_, SEXP gamma_, SEXP group_multiplier, SEXP dfmax_, SEXP gmax_, SEXP warn_, SEXP user_) {

// Lengths/dimensions
int n = length(d_);
Expand All @@ -55,6 +55,9 @@ SEXP gdfit_cox(SEXP X_, SEXP d_, SEXP penalty_, SEXP K1_, SEXP K0_, SEXP lambda,

// Pointers
double *X = REAL(X_);
double *start_time = REAL(start_);
double *start_o = REAL(start_o_);
double *stop_time = REAL(stop_);
double *d = REAL(d_);
const char *penalty = CHAR(STRING_ELT(penalty_, 0));
int *K1 = INTEGER(K1_);
Expand Down Expand Up @@ -98,6 +101,8 @@ SEXP gdfit_cox(SEXP X_, SEXP d_, SEXP penalty_, SEXP K1_, SEXP K0_, SEXP lambda,
for (int g=0; g<J; g++) e[g] = 0;
int lstart, ng, nv, violations;
double shift, l1, l2, nullDev, v, s, maxChange;
double current_sum;
int start_idx, stop_idx, ss;

// Initialization
// If lam[0]=lam_max, skip lam[0] -- closed form sol'n available
Expand Down Expand Up @@ -147,6 +152,19 @@ SEXP gdfit_cox(SEXP X_, SEXP d_, SEXP penalty_, SEXP K1_, SEXP K0_, SEXP lambda,
for (int i=n-2; i>=0; i--) {
rsk[i] = rsk[i+1] + haz[i];
}
current_sum = 0;
start_idx = n;
stop_idx = n;
while(start_idx>0 && stop_idx>0){
ss = start_o[start_idx];
if(start_time[ss] < stop_time[stop_idx]){
rsk[stop_idx] = rsk[stop_idx] - current_sum;
stop_idx = stop_idx - 1;
}else{
current_sum = current_sum + haz[ss];
start_idx = start_idx - 1;
}
}
for (int i=0; i<n; i++) {
REAL(Loss)[l] += d[i]*eta[i] - d[i]*log(rsk[i]);
}
Expand Down
8 changes: 4 additions & 4 deletions src/grpreg_init.c
Original file line number Diff line number Diff line change
Expand Up @@ -7,21 +7,21 @@

/* .Call calls */
extern SEXP gdfit_glm(SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP);
extern SEXP gdfit_cox(SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP);
extern SEXP gdfit_cox(SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP);
extern SEXP gdfit_gaussian(SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP);
extern SEXP lcdfit_glm(SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP);
extern SEXP lcdfit_cox(SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP);
extern SEXP lcdfit_cox(SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP);
extern SEXP lcdfit_gaussian(SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP, SEXP);
extern SEXP maxgrad(SEXP, SEXP, SEXP, SEXP);
extern SEXP maxprod(SEXP, SEXP, SEXP, SEXP);
extern SEXP standardize(SEXP);

static const R_CallMethodDef CallEntries[] = {
{"gdfit_glm", (DL_FUNC) &gdfit_glm, 16},
{"gdfit_cox", (DL_FUNC) &gdfit_cox, 15},
{"gdfit_cox", (DL_FUNC) &gdfit_cox, 18},
{"gdfit_gaussian", (DL_FUNC) &gdfit_gaussian, 15},
{"lcdfit_glm", (DL_FUNC) &lcdfit_glm, 18},
{"lcdfit_cox", (DL_FUNC) &lcdfit_cox, 17},
{"lcdfit_cox", (DL_FUNC) &lcdfit_cox, 20},
{"lcdfit_gaussian", (DL_FUNC) &lcdfit_gaussian, 16},
{"maxgrad", (DL_FUNC) &maxgrad, 4},
{"maxprod", (DL_FUNC) &maxprod, 4},
Expand Down
22 changes: 20 additions & 2 deletions src/lcdfit_cox.c
Original file line number Diff line number Diff line change
Expand Up @@ -133,7 +133,7 @@ int gLCD_cCheck(double *b, const char *penalty, double *X, double *r, double *et
return(violations);
}

SEXP lcdfit_cox(SEXP X_, SEXP d_, SEXP penalty_, SEXP K1_, SEXP K0_,
SEXP lcdfit_cox(SEXP X_, SEXP start_, SEXP start_o_, SEXP stop_, SEXP d_, SEXP penalty_, SEXP K1_, SEXP K0_,
SEXP lambda, SEXP alpha_, SEXP eps_, SEXP delta_, SEXP gamma_,
SEXP tau_, SEXP max_iter_, SEXP group_multiplier, SEXP dfmax_,
SEXP gmax_, SEXP warn_, SEXP user_) {
Expand All @@ -146,6 +146,9 @@ SEXP lcdfit_cox(SEXP X_, SEXP d_, SEXP penalty_, SEXP K1_, SEXP K0_,

// Pointers
double *X = REAL(X_);
double *start_time = REAL(start_);
double *start_o = REAL(start_o_);
double *stop_time = REAL(stop_);
double *d = REAL(d_);
const char *penalty = CHAR(STRING_ELT(penalty_, 0));
int *K1 = INTEGER(K1_);
Expand Down Expand Up @@ -191,7 +194,9 @@ SEXP lcdfit_cox(SEXP X_, SEXP d_, SEXP penalty_, SEXP K1_, SEXP K0_,
for (int j=0; j<p; j++) e[j] = 0;
int lstart, ng, nv, violations;
double shift, l1, l2, nullDev, u, v, s, xwr, xwx, maxChange;

double current_sum;
int start_idx, stop_idx, ss;

// Initialization
rsk[n-1] = 1;
for (int i=n-2; i>=0; i--) rsk[i] = rsk[i+1] + 1;
Expand Down Expand Up @@ -256,6 +261,19 @@ SEXP lcdfit_cox(SEXP X_, SEXP d_, SEXP penalty_, SEXP K1_, SEXP K0_,
for (int i=n-2; i>=0; i--) {
rsk[i] = rsk[i+1] + haz[i];
}
current_sum = 0;
start_idx = n;
stop_idx = n;
while(start_idx>0 && stop_idx>0){
ss = start_o[start_idx];
if(start_time[ss] < stop_time[stop_idx]){
rsk[stop_idx] = rsk[stop_idx] - current_sum;
stop_idx = stop_idx - 1;
}else{
current_sum = current_sum + haz[ss];
start_idx = start_idx - 1;
}
}
for (int i=0; i<n; i++) {
REAL(Loss)[l] += d[i]*eta[i] - d[i]*log(rsk[i]);
}
Expand Down