Title: | Fine-Gray Regression via Forward-Backward Scan |
---|---|
Description: | In competing risks regression, the proportional subdistribution hazards (PSH) model is popular for its direct assessment of covariate effects on the cumulative incidence function. This package allows for both penalized and unpenalized PSH regression in linear time using a novel forward-backward scan. Penalties include Ridge, Lease Absolute Shrinkage and Selection Operator (LASSO), Smoothly Clipped Absolute Deviation (SCAD), Minimax Concave Plus (MCP), and elastic net <doi: 10.32614/RJ-2021-010>. |
Authors: | Eric S. Kawaguchi [aut, cre] |
Maintainer: | Eric S. Kawaguchi <[email protected]> |
License: | GPL-3 |
Version: | 1.24.10 |
Built: | 2024-11-02 05:23:29 UTC |
Source: | https://github.com/erickawaguchi/fastcmprsk |
Similar functional utility to coef
methods.
## S3 method for class 'fcrr' AIC(object, ..., k = 2)
## S3 method for class 'fcrr' AIC(object, ..., k = 2)
object |
|
... |
Additional arguments. Not implemented. |
k |
Numeric, the penalty per parameter to be used; the default k = 2 is the classical AIC. |
A numeric value with the corresponding AIC (or BIC, or ..., depending on k
).
Similar functional utility to coef
methods.
## S3 method for class 'fcrrp' AIC(object, ..., k = 2)
## S3 method for class 'fcrrp' AIC(object, ..., k = 2)
object |
|
... |
Additional arguments. Not implemented. |
k |
Numeric, the penalty per parameter to be used; the default k = 2 is the classical AIC. |
A numeric value with the corresponding AIC (or BIC, or ..., depending on k
).
Similar functional utility to coef
methods.
## S3 method for class 'fcrr' coef(object, ...)
## S3 method for class 'fcrr' coef(object, ...)
object |
|
... |
Additional arguments. Not implemented. |
Coefficients extracted from the model object object
.
Similar functional utility to coef
methods.
## S3 method for class 'fcrrp' coef(object, ...)
## S3 method for class 'fcrrp' coef(object, ...)
object |
|
... |
Additional arguments. Not implemented. |
Coefficients extracted from the model object object
.
Computes confidence intervals for one or more parameters in a fitted model of class fcrr
.
## S3 method for class 'fcrr' confint(object, parm, level = 0.95, digits = max(options()$digits - 5, 2), ...)
## S3 method for class 'fcrr' confint(object, parm, level = 0.95, digits = max(options()$digits - 5, 2), ...)
object |
|
parm |
a specification of which parameters are to be given confidence intervals, either a vector of numbers or a vector of names. If missing, all parameters are considered. |
level |
the confidence level required |
digits |
Number of significant difits to round to. |
... |
Additional arguments. Not implemented. |
Prints out table of confidence intervals for the Fine-Gray model.
A matrix (or vector) with columns giving lower and upper confidence limits for each coefficient estimate.
Create a competing risk object, used as a response variable in the model formula for fastCrr
and fastCrrp
.
Adapted from the Surv
object.
Crisk(ftime, fstatus, cencode = 0, failcode = 1, silent = TRUE)
Crisk(ftime, fstatus, cencode = 0, failcode = 1, silent = TRUE)
ftime |
A vector of event/censoring times. |
fstatus |
A vector with unique code for each event type and a separate code for censored observations. |
cencode |
Integer: code of |
failcode |
Integer: code of |
silent |
Logical: print information about coding. |
Returns an object, used as a response variable, of class Crisk
.
time |
vector of observed event times |
status |
vector of event indicators. 0 = censored, 1 = event of interest, 2 = competing risks |
Fine J. and Gray R. (1999) A proportional hazards model for the subdistribution of a competing risk. JASA 94:496-509.
Surv
library(fastcmprsk) set.seed(10) ftime <- rexp(200) fstatus <- sample(0:2, 200, replace = TRUE) obj <- Crisk(ftime, fstatus, silent = FALSE)
library(fastcmprsk) set.seed(10) ftime <- rexp(200) fstatus <- sample(0:2, 200, replace = TRUE) obj <- Crisk(ftime, fstatus, silent = FALSE)
Estimates parameters for the proportional subdistribution hazards model using two-way linear scan approach.
fastCrr( formula, data, eps = 1e-06, max.iter = 1000, getBreslowJumps = TRUE, standardize = TRUE, variance = TRUE, var.control = varianceControl(B = 100, useMultipleCores = FALSE), returnDataFrame = FALSE )
fastCrr( formula, data, eps = 1e-06, max.iter = 1000, getBreslowJumps = TRUE, standardize = TRUE, variance = TRUE, var.control = varianceControl(B = 100, useMultipleCores = FALSE), returnDataFrame = FALSE )
formula |
a formula object, with the response on the left of a ~ operator, and the terms on the right. The response must be a Crisk object as returned by the |
data |
a data.frame in which to interpret the variables named in the formula. |
eps |
Numeric: algorithm stops when the relative change in any coefficient is less than |
max.iter |
Numeric: maximum iterations to achieve convergence (default is 1000) |
getBreslowJumps |
Logical: Output jumps in Breslow estimator for the cumulative hazard. |
standardize |
Logical: Standardize design matrix. |
variance |
Logical: Get standard error estimates for parameter estimates via bootstrap. |
var.control |
List of options for variance estimation. |
returnDataFrame |
Logical: Return (ordered) data frame. |
Fits the 'proportional subdistribution hazards' regression model described in Fine and Gray (1999) using a novel two-way linear scan approach.
By default, the Crisk
object will specify which observations are censored (0), the event of interest (1), or competing risks (2).
Returns a list of class fcrr
.
coef |
the estimated regression coefficients |
var |
estimated variance-covariance matrix via bootstrap (if |
logLik |
log-pseudo likelihood at the estimated regression coefficients |
logLik.null |
log-pseudo likelihood when the regression coefficients are 0 |
lrt |
log-pseudo likelihood ratio test statistic for the estimated model vs. the null model. |
iter |
iterations of coordinate descent until convergence |
converged |
logical. |
breslowJump |
Jumps in the Breslow baseline cumulative hazard (used by |
uftime |
vector of unique failure (event) times |
isVariance |
logical to return if variance is chosen to be estimated |
df |
returned ordered data frame if |
Fine J. and Gray R. (1999) A proportional hazards model for the subdistribution of a competing risk. JASA 94:496-509.
#' Kawaguchi, E.S., Shen J.I., Suchard, M. A., Li, G. (2020) Scalable Algorithms for Large Competing Risks Data, Journal of Computational and Graphical Statistics
library(fastcmprsk) set.seed(10) ftime <- rexp(200) fstatus <- sample(0:2, 200, replace = TRUE) cov <- matrix(runif(1000), nrow = 200) dimnames(cov)[[2]] <- c('x1','x2','x3','x4','x5') fit <- fastCrr(Crisk(ftime, fstatus) ~ cov, variance = FALSE) # Not run: How to set up multiple cores for boostrapping # library(doParallel) # make sure necessary packages are loaded # myClust <- makeCluster(2) # registerDoParallel(myClust) # fit1 <- fastCrr(Crisk(ftime, fstatus) ~ cov, variance = TRUE, # var.control = varianceControl(B = 100, useMultipleCores = TRUE)) # stopCluster(myClust)
library(fastcmprsk) set.seed(10) ftime <- rexp(200) fstatus <- sample(0:2, 200, replace = TRUE) cov <- matrix(runif(1000), nrow = 200) dimnames(cov)[[2]] <- c('x1','x2','x3','x4','x5') fit <- fastCrr(Crisk(ftime, fstatus) ~ cov, variance = FALSE) # Not run: How to set up multiple cores for boostrapping # library(doParallel) # make sure necessary packages are loaded # myClust <- makeCluster(2) # registerDoParallel(myClust) # fit1 <- fastCrr(Crisk(ftime, fstatus) ~ cov, variance = TRUE, # var.control = varianceControl(B = 100, useMultipleCores = TRUE)) # stopCluster(myClust)
Performs penalized regression for the proportional subdistribution hazards model. Penalties currently include LASSO, MCP, SCAD, and ridge regression. User-specificed weights can be assigned to the penalty for each coefficient (e.g. implementing adaptive LASSO and broken adaptive ridge regerssion).
fastCrrp( formula, data, eps = 1e-06, max.iter = 1000, getBreslowJumps = TRUE, standardize = TRUE, penalty = c("LASSO", "RIDGE", "MCP", "SCAD", "ENET"), lambda = NULL, alpha = 0, lambda.min.ratio = 0.001, nlambda = 25, penalty.factor, gamma = switch(penalty, scad = 3.7, 2.7) )
fastCrrp( formula, data, eps = 1e-06, max.iter = 1000, getBreslowJumps = TRUE, standardize = TRUE, penalty = c("LASSO", "RIDGE", "MCP", "SCAD", "ENET"), lambda = NULL, alpha = 0, lambda.min.ratio = 0.001, nlambda = 25, penalty.factor, gamma = switch(penalty, scad = 3.7, 2.7) )
formula |
a formula object, with the response on the left of a ~ operator, and the terms on the right. The response must be a Crisk object as returned by the |
data |
a data.frame in which to interpret the variables named in the formula. |
eps |
Numeric: algorithm stops when the relative change in any coefficient is less than |
max.iter |
Numeric: maximum iterations to achieve convergence (default is 1000) |
getBreslowJumps |
Logical: Output jumps in Breslow estimator for the cumulative hazard (one for each value of lambda). |
standardize |
Logical: Standardize design matrix. |
penalty |
Character: Penalty to be applied to the model. Options are "lasso", "scad", "ridge", "mcp", and "enet". |
lambda |
A user-specified sequence of |
alpha |
L1/L2 weight for elastic net regression. |
lambda.min.ratio |
Smallest value for |
nlambda |
Number of |
penalty.factor |
A vector of weights applied to the penalty for each coefficient. Vector must be of length equal to the number of columns in |
gamma |
Tuning parameter for the MCP/SCAD penalty. Default is 2.7 for MCP and 3.7 for SCAD and should be left unchanged. |
The fastCrrp
functions performed penalized Fine-Gray regression.
Parameter estimation is performed via cyclic coordinate descent and using a two-way linear scan approach to efficiently
calculate the gradient and Hessian values. Current implementation includes LASSO, SCAD, MCP, and ridge regression.
Returns a list of class fcrrp
.
coef |
fitted coefficients matrix with |
logLik |
vector of log-pseudo likelihood at the estimated regression coefficients |
logLik.null |
log-pseudo likelihood when the regression coefficients are 0 |
lambda.path |
sequence of tuning parameter values |
iter |
number of iterations needed until convergence at each tuning parameter value |
converged |
convergence status at each tuning parameter value |
breslowJump |
Jumps in the Breslow baseline cumulative hazard (used by |
uftime |
vector of unique failure (event) times |
penalty |
same as above |
gamma |
same as above |
above |
same as above |
Fu, Z., Parikh, C.R., Zhou, B. (2017) Penalized variable selection in competing risks regression. Lifetime Data Analysis 23:353-376.
Breheny, P. and Huang, J. (2011) Coordinate descent algorithms for nonconvex penalized regression, with applications to biological feature selection. Ann. Appl. Statist., 5: 232-253.
Fine J. and Gray R. (1999) A proportional hazards model for the subdistribution of a competing risk. JASA 94:496-509.
Kawaguchi, E.S., Shen J.I., Suchard, M. A., Li, G. (2020) Scalable Algorithms for Large Competing Risks Data, Journal of Computational and Graphical Statistics
library(fastcmprsk) set.seed(10) ftime <- rexp(200) fstatus <- sample(0:2, 200, replace = TRUE) cov <- matrix(runif(1000), nrow = 200) dimnames(cov)[[2]] <- c('x1','x2','x3','x4','x5') fit <- fastCrrp(Crisk(ftime, fstatus) ~ cov, lambda = 1, penalty = "RIDGE") fit$coef
library(fastcmprsk) set.seed(10) ftime <- rexp(200) fstatus <- sample(0:2, 200, replace = TRUE) cov <- matrix(runif(1000), nrow = 200) dimnames(cov)[[2]] <- c('x1','x2','x3','x4','x5') fit <- fastCrrp(Crisk(ftime, fstatus) ~ cov, lambda = 1, penalty = "RIDGE") fit$coef
Similar functional utility to coef
methods.
## S3 method for class 'fcrr' logLik(object, ...)
## S3 method for class 'fcrr' logLik(object, ...)
object |
|
... |
Additional arguments. Not implemented. |
Returns the log-pseudo likelihood of object object
.
Similar functional utility to coef
methods.
## S3 method for class 'fcrrp' logLik(object, ...)
## S3 method for class 'fcrrp' logLik(object, ...)
object |
|
... |
Additional arguments. Not implemented. |
Returns the log-pseudo likelihood of object object
.
Plots solution path for penalized methods
## S3 method for class 'fcrrp' plot(x, ...)
## S3 method for class 'fcrrp' plot(x, ...)
x |
|
... |
additional arguments to |
Plots solution path for penalized methods. x-axis: log tuning parameter values. y-axis: coeffcient estimates.
A plot of the solution path for the chosen penalized method.
Plots predicted cumulative incidence function
## S3 method for class 'predict.fcrr' plot(x, ...)
## S3 method for class 'predict.fcrr' plot(x, ...)
x |
|
... |
additional arguments to |
A plot of the estimated cumulative incidence function.
Predicts cumulative incidence function from a fcrr
object.
## S3 method for class 'fcrr' predict( object, newdata, getBootstrapVariance = TRUE, var.control = varianceControl(B = 100, useMultipleCores = FALSE), type = "none", alpha = 0.05, tL = NULL, tU = NULL, ... )
## S3 method for class 'fcrr' predict( object, newdata, getBootstrapVariance = TRUE, var.control = varianceControl(B = 100, useMultipleCores = FALSE), type = "none", alpha = 0.05, tL = NULL, tU = NULL, ... )
object |
Output from |
newdata |
A set of covariate values to predict the CIF. |
getBootstrapVariance |
Logical: Calculate variance for CIF via bootstrap. |
var.control |
List of variance parameters from |
type |
Confidence intervals or confidence bands. |
alpha |
Significance level to compute intervals or bands |
tL |
Lower time for band estimation. |
tU |
Upper time for band estimation. |
... |
additional arguments affecting the fastCrr procedure. |
B |
Number of bootstrap samples for variance estimation. |
Calculates the CIF using fcrr
output conditional on newdata
.
Returns a list of class predict.fcrr
.
ftime |
Unique observed failure times |
CIF |
predicted CIF at time |
lower |
lower interval/band limit |
upper |
upper interval/band limit |
type |
same as original argument |
Fine J. and Gray R. (1999) A proportional hazards model for the subdistribution of a competing risk. JASA 94:496-509.
library(fastcmprsk) set.seed(10) ftime <- rexp(200) fstatus <- sample(0:2, 200, replace = TRUE) cov <- matrix(runif(1000), nrow = 200) dimnames(cov)[[2]] <- c('x1','x2','x3','x4','x5') fit <- fastCrr(Crisk(ftime, fstatus) ~ cov, returnDataFrame = TRUE) cov2 <- rnorm(5) predict(fit, newdata = cov2)
library(fastcmprsk) set.seed(10) ftime <- rexp(200) fstatus <- sample(0:2, 200, replace = TRUE) cov <- matrix(runif(1000), nrow = 200) dimnames(cov)[[2]] <- c('x1','x2','x3','x4','x5') fit <- fastCrr(Crisk(ftime, fstatus) ~ cov, returnDataFrame = TRUE) cov2 <- rnorm(5) predict(fit, newdata = cov2)
Prints summary statistics of a fcrr x
## S3 method for class 'summary.fcrr' print(x, digits = max(options()$digits - 4, 3), ...)
## S3 method for class 'summary.fcrr' print(x, digits = max(options()$digits - 4, 3), ...)
x |
output from |
digits |
digits for rounding. |
... |
additional arguments to |
Prints the convergence status, log-pseudo likelihood, the estimated coefficients, the estimated standard errors, and the two-sided p-values for the test of the individual coefficients equal to 0.
Prints the convergence status, log-pseudo likelihood, the estimated coefficients, the estimated standard errors, and the two-sided p-values for the test of the individual coefficients equal to 0.
Simulate data from the model proposed in Fine and Gray (1999) for two causes. Cause 1 is assumed to be of primary importance.
simulateTwoCauseFineGrayModel( nobs, beta1, beta2, X = NULL, u.min = 0, u.max, p = 0.5, returnX = FALSE )
simulateTwoCauseFineGrayModel( nobs, beta1, beta2, X = NULL, u.min = 0, u.max, p = 0.5, returnX = FALSE )
nobs |
Integer: Number of observations in simulated dataset. |
beta1 |
A vector of effect sizes for cause 1 of length ncovs |
beta2 |
A vector of effect sizes for cause 2 of length ncovs |
X |
A matrix of fixed covariates (nobs x ncovs). If |
u.min |
Numeric: controls lower bound of censoring distribution where C ~ U(u.min, u.max) |
u.max |
Numeric: controls upper bound of censoring distribution where C ~ U(u.min, u.max) |
p |
Numeric: value between 0 and 1 which controls the mixture probability. |
returnX |
Logical: Whether to return |
The function simulates data according to the setup by Fine and Gray (1999). See their paper for more information.
Returns a list with the following:
ftime |
vector of |
ftime |
vector of |
X |
design matrix if |
Fine J. and Gray R. (1999) A proportional hazards model for the subdistribution of a competing risk. JASA 94:496-509.
set.seed(2019) nobs <- 500 beta1 <- c(0.40, -0.40, 0, -0.50, 0, 0.60, 0.75, 0, 0, -0.80) beta2 <- -beta1 Z <- matrix(rnorm(nobs * length(beta1)), nrow = nobs) dat <- simulateTwoCauseFineGrayModel(nobs, beta1, beta2, Z, u.min = 0, u.max = 1, p = 0.5)
set.seed(2019) nobs <- 500 beta1 <- c(0.40, -0.40, 0, -0.50, 0, 0.60, 0.75, 0, 0, -0.80) beta2 <- -beta1 Z <- matrix(rnorm(nobs * length(beta1)), nrow = nobs) dat <- simulateTwoCauseFineGrayModel(nobs, beta1, beta2, Z, u.min = 0, u.max = 1, p = 0.5)
Generate and print summaries of fastCrr
output.
## S3 method for class 'fcrr' summary( object, conf.int = TRUE, alpha = 0.05, digits = max(options()$digits - 5, 2), ... )
## S3 method for class 'fcrr' summary( object, conf.int = TRUE, alpha = 0.05, digits = max(options()$digits - 5, 2), ... )
object |
|
conf.int |
Logical. Whether or not to outut confidence intervals. |
alpha |
Significance level of the confidence intervals. |
digits |
Numer of significant difits to round to. |
... |
additional arguments to |
The summary method produces an ANOVA table for the coefficient estimates of the Fine-Gray model.
The form of the value returned by summary
depends on the class of its argument. See the documentation of the particular methods for details of what is produced by that method.
Controls for variance calculation for the fastcmprsk package.
varianceControl( B = 100L, seed = 1991L, useMultipleCores = FALSE, extractMatrix = FALSE )
varianceControl( B = 100L, seed = 1991L, useMultipleCores = FALSE, extractMatrix = FALSE )
B |
Integer: Number of bootstrap samples needed for variance estimation. |
seed |
Integer: Seed value for bootstrapping. Results may differ if parallelized. |
useMultipleCores |
Logical: Set to TRUE if parallelizing. (Default is FALSE). |
extractMatrix |
Logical: Extract matrix of bootstrap estimates (Default is FALSE) |
Variance-covariance estimation is done via bootstrap.
Independent bootstrap runs can be performed both in serial and parallel. Parallelization is done via the
doParallel
package.
Returns a list for variance options inputted into fastCrr
.
B |
same as what is defined in argument. |
seed |
same as what is defined in argument. |
mcores |
same as what is defined in argument |
extract |
same as what is defined in argument |
library(fastcmprsk) set.seed(10) ftime <- rexp(200) fstatus <- sample(0:2, 200, replace = TRUE) cov <- matrix(runif(1000), nrow = 200) dimnames(cov)[[2]] <- c('x1','x2','x3','x4','x5') vc <- varianceControl(B = 100, seed = 2019, useMultipleCores = FALSE) fit1 <- fastCrr(Crisk(ftime, fstatus) ~ cov, variance = TRUE, var.control = vc) fit1$var # Estimated covariance matrix via bootstrap
library(fastcmprsk) set.seed(10) ftime <- rexp(200) fstatus <- sample(0:2, 200, replace = TRUE) cov <- matrix(runif(1000), nrow = 200) dimnames(cov)[[2]] <- c('x1','x2','x3','x4','x5') vc <- varianceControl(B = 100, seed = 2019, useMultipleCores = FALSE) fit1 <- fastCrr(Crisk(ftime, fstatus) ~ cov, variance = TRUE, var.control = vc) fit1$var # Estimated covariance matrix via bootstrap
Similar functional utility to vcov
methods.
## S3 method for class 'fcrr' vcov(object, ...)
## S3 method for class 'fcrr' vcov(object, ...)
object |
|
... |
Additional arguments. Not implemented. |
Returns the estimated variance-covariance matrix (via bootstrap) from object object
.