Title: | Weibull Regression for a Right-Censored Endpoint with Interval-Censored Covariate |
---|---|
Description: | The function SurvRegCens() of this package allows estimation of a Weibull Regression for a right-censored endpoint, one interval-censored covariate, and an arbitrary number of non-censored covariates. Additional functions allow to switch between different parametrizations of Weibull regression used by different R functions, inference for the mean difference of two arbitrarily censored Normal samples, and estimation of canonical parameters from censored samples for several distributional assumptions. Hubeaux, S. and Rufibach, K. (2014) <arXiv:1402.0432>. |
Authors: | Stanislas Hubeaux <[email protected]> and Kaspar Rufibach <[email protected]> |
Maintainer: | Stanislas Hubeaux <[email protected]> |
License: | GPL (>= 2) |
Version: | 1.7 |
Built: | 2025-01-19 02:58:58 UTC |
Source: | https://github.com/cran/SurvRegCensCov |
The function SurvRegCens
of this package allows estimation of a Weibull Regression for a right-censored endpoint, one interval-censored covariate, and an arbitrary number of non-censored covariates. Additional functions allow to switch between different parametrizations of Weibull regression used by different R
functions (ConvertWeibull
, WeibullReg
, WeibullDiag
), inference for the mean difference of two arbitrarily censored Normal samples (NormalMeanDiffCens
), and estimation of canonical parameters from censored samples for several distributional assumptions (ParamSampleCens
).
Package: | SurvRegCensCov |
Type: | Package |
Version: | 1.7 |
Date: | 2023-09-27 |
License: | GPL (>=2) |
Stanislas Hubeaux (maintainer), [email protected]
Kaspar Rufibach, [email protected]
http://www.kasparrufibach.ch
We thank Sarah Haile for contributing the functions ConvertWeibull
, WeibullReg
, WeibullDiag
to the package.
Hubeaux, S. (2013). Estimation from left- and/or interval-censored samples. Technical report, Biostatistics Oncology, F. Hoffmann-La Roche Ltd.
Hubeaux, S. (2013). Parametric Surival Regression Model with left- and/or interval-censored covariate. Technical report, Biostatistics Oncology, F. Hoffmann-La Roche Ltd.
Hubeaux, S. and Rufibach, K. (2014). SurvRegCensCov: Weibull Regression for a Right-Censored Endpoint with a Censored Covariate. Preprint, https://arxiv.org/abs/1402.0432.
Lynn, H. S. (2001). Maximum likelihood inference for left-censored HIV RNA data. Stat. Med., 20, 33–45.
Sattar, A., Sinha, S. K. and Morris, N. J. (2012). A Parametric Survival Model When a Covariate is Subject to Left-Censoring. Biometrics & Biostatistics, S3(2).
# The main functions in this package are illustrated in their respective help files.
# The main functions in this package are illustrated in their respective help files.
Evaluates the cumulative distribution function using the integral of its density function.
CDF(c, density)
CDF(c, density)
c |
Value at which the CDF is to be evaluated. |
density |
Density function. |
Function not intended to be invoked by the user.
Stanislas Hubeaux, [email protected]
Kaspar Rufibach, [email protected]
http://www.kasparrufibach.ch
Given a vector of realizations of a continuous random variable, interval-, left-, or right-censor these numbers at given boundaries. Useful when setting up simulations involving censored observations.
censorContVar(x, LLOD = NA, ULOD = NA)
censorContVar(x, LLOD = NA, ULOD = NA)
x |
Vector of random numbers. |
LLOD |
Lower limit where |
ULOD |
Upper limit where |
A data.frame
as specified by code = interval2
in Surv
.
Stanislas Hubeaux, [email protected]
Kaspar Rufibach, [email protected]
http://www.kasparrufibach.ch
## random vector x <- rnorm(200) ## interval-censor this vector at -1 and 0.5 censorContVar(x, -1, 0.5)
## random vector x <- rnorm(200) ## interval-censor this vector at -1 and 0.5 censorContVar(x, -1, 0.5)
coef
method for class "src"
.
## S3 method for class 'src' coef(object, ...)
## S3 method for class 'src' coef(object, ...)
object |
An object of class |
... |
Further arguments. |
The function coef.src
returns the estimated parameters of the Weibull regression when calling SurvRegCens
.
Stanislas Hubeaux, [email protected]
Kaspar Rufibach, [email protected]
http://www.kasparrufibach.ch
Hubeaux, S. (2013). Parametric Surival Regression Model with left- and/or interval-censored covariate. Technical report, Biostatistics Oncology, F. Hoffmann-La Roche Ltd.
Hubeaux, S. and Rufibach, K. (2014). SurvRegCensCov: Weibull Regression for a Right-Censored Endpoint with a Censored Covariate. Preprint, https://arxiv.org/abs/1402.0432.
Sattar, A., Sinha, S. K. and Morris, N. J. (2012). A Parametric Survival Model When a Covariate is Subject to Left-Censoring. Biometrics & Biostatistics, S3(2).
## See help file of function "SurvRegCens".
## See help file of function "SurvRegCens".
Transforms output from survreg
using the Weibull distribution to a more natural parameterization. See details and the vignette for more information.
ConvertWeibull(model, conf.level = 0.95)
ConvertWeibull(model, conf.level = 0.95)
model |
A |
conf.level |
Confidence level used to produce two-sided |
The survreg
function fits a Weibull accelerated failure time model of the form
where is a matrix of covariates, and
has the extreme value distribution,
is the intercept,
is a vector of parameters for each of the covariates, and
is the scale. The usual
parameterization of the model, however, is defined by hazard function
The transformation is as follows: ,
, and
, and estimates of the standard errors can be found using the delta method.
The Weibull distribution has the advantage of having two separate interpretations. The first, via proportional hazards, leads to a hazard ratio, defined by . The second, of accelerated failure times, leads to an event time ratio (also known as an acceleration factor), defined by
.
Further details regarding the transformations of the parameters and their standard errors can be found in Klein and
Moeschberger (2003, Chapter 12). An explanation of event time ratios for the accelerated failure time interpretation of the model can be found in Carroll (2003). A general overview can be found in the vignette("weibull")
of this package.
vars |
A matrix containing the values of the transformed parameters and their standard errors |
HR |
A matrix containing the hazard ratios for the covariates, and |
ETR |
A matrix containing the event time ratios for the covariates, and |
Sarah R. Haile, Epidemiology, Biostatistics and Prevention Institute (EBPI), University of Zurich, [email protected]
Carroll, K. (2003). On the use and utility of the Weibull model in the analysis of survival data. Controlled Clinical Trials, 24, 682–701.
Klein, J. and Moeschberger, M. (2003). Survival analysis: techniques for censored and truncated data. 2nd edition, Springer.
This function is used by WeibullReg
.
data(larynx) ConvertWeibull(survreg(Surv(time, death) ~ stage + age, larynx), conf.level = 0.95)
data(larynx) ConvertWeibull(survreg(Surv(time, death) ~ stage + age, larynx), conf.level = 0.95)
A study of 90 males with laryngeal cancer was performed, comparing survival times. Each patient's age, year of diagnosis, and disease stage was noted, see Kardaun (1983) and Klein and Moeschberger (2003).
data(larynx)
data(larynx)
A data frame with 90 observations on the following 5 variables.
stage
Disease stage (1-4) from TNM cancer staging classification.
time
Time from first treatment until death, or end of study.
age
Age at diagnosis.
year
Year of diagnosis.
death
Indicator of death [1, if patient died at time t; 0, otherwise].
https://www.mcw.edu/-/media/MCW/Departments/Biostatistics/datafromsection18.txt?la=en
Kardaun, O. (1983). Statistical survival analysis of male larynx-cancer patients-a case study. Statistica Neerlandica, 37, 103–125.
Klein, J. and Moeschberger, M. (2003). Survival analysis: techniques for censored and truncated data. 2nd edition, Springer.
library(survival) data(larynx) Surv(larynx$time, larynx$death)
library(survival) data(larynx) Surv(larynx$time, larynx$death)
logLik
method for class "src"
.
## S3 method for class 'src' logLik(object, ...)
## S3 method for class 'src' logLik(object, ...)
object |
An object of class |
... |
Further arguments. |
The function logLik.src
returns the value of the log-likelihood at the maximum likelihood estimate, as well as the corresponding degrees of freedom.
Stanislas Hubeaux, [email protected]
Kaspar Rufibach, [email protected]
http://www.kasparrufibach.ch
Hubeaux, S. (2013). Parametric Surival Regression Model with left- and/or interval-censored covariate. Technical report, Biostatistics Oncology, F. Hoffmann-La Roche Ltd.
Hubeaux, S. and Rufibach, K. (2014). SurvRegCensCov: Weibull Regression for a Right-Censored Endpoint with a Censored Covariate. Preprint, https://arxiv.org/abs/1402.0432.
Sattar, A., Sinha, S. K. and Morris, N. J. (2012). A Parametric Survival Model When a Covariate is Subject to Left-Censoring. Biometrics & Biostatistics, S3(2).
## See help file of function "SurvRegCens".
## See help file of function "SurvRegCens".
Computes the log-likelihood function for a censored sample, according to a specified distributional assumptions. Available distributions are Normal, Weibull, Logistic, and Gamma.
LoglikNormalCens(x, data, lowerbound, vdelta) LoglikWeibullCens(x, data, lowerbound, vdelta) LoglikLogisticCens(x, data, lowerbound, vdelta) LoglikGammaCens(x, data, lowerbound, vdelta)
LoglikNormalCens(x, data, lowerbound, vdelta) LoglikWeibullCens(x, data, lowerbound, vdelta) LoglikLogisticCens(x, data, lowerbound, vdelta) LoglikGammaCens(x, data, lowerbound, vdelta)
x |
Two-dimensional vector giving the canonical parameters of the distribution. |
data |
Observed or censored event times. |
lowerbound |
A vector that collect lower bounds for the interval-censored observations. If no lower bound is available then put |
vdelta |
A vector which indicates censoring (0: censored, 1: not censored). |
Function not intended to be invoked by the user.
Stanislas Hubeaux, [email protected]
Kaspar Rufibach, [email protected]
http://www.kasparrufibach.ch
Hubeaux, S. (2013). Estimation from left- and/or interval-censored samples. Technical report, Biostatistics Oncology, F. Hoffmann-La Roche Ltd.
Lynn, H. S. (2001). Maximum likelihood inference for left-censored HIV RNA data. Stat. Med., 20, 33–45.
Reparametrization of the log likelihood function for a normally distributed censored sample such that the mean difference is a parameter of the function, thus allowing to be made inference on. The mean difference is computed as sample 1 - sample 2.
LoglikNormalDeltaCens(x, data1, lowerbound1, vdelta1, data2, lowerbound2, vdelta2)
LoglikNormalDeltaCens(x, data1, lowerbound1, vdelta1, data2, lowerbound2, vdelta2)
x |
A vector of four components where the first component corresponds to the mean of the normal distribution of the first sample, the second component corresponds to mean difference between the two samples: sample 1 - sample 2, the third component corresponds to the standard deviation of the normal distribution of the first sample, and the fourth component corresponds to the standard deviation of the normal distribution of the second sample. |
data1 |
A vector of data corresponding to the first sample. |
lowerbound1 |
A vector which corresponds to the lower bounds for the interval-censored observations of the vector of data corresponding to the first sample. If no lower bound is available then put |
vdelta1 |
A vector which indicates for censoring for the first sample (0: censored, 1: not censored). |
data2 |
A vector of data corresponding to the second sample. |
lowerbound2 |
A vector which corresponds to the lower bounds for the interval-censored observations of the vector of data corresponding to the second sample. If no lower bound is available then put |
vdelta2 |
A vector which indicates for censoring for the second sample (0: censored, 1: not censored). |
Function not intended to be invoked by the user.
Stanislas Hubeaux, [email protected]
Kaspar Rufibach, [email protected]
http://www.kasparrufibach.ch
Hubeaux, S. (2013). Estimation from left- and/or interval-censored samples. Technical report, Biostatistics Oncology, F. Hoffmann-La Roche Ltd.
Lynn, H. S. (2001). Maximum likelihood inference for left-censored HIV RNA data. Stat. Med., 20, 33–45.
Computes the log-likelihood function of a Weibull Survival Regression Model allowing for an interval-censored covariate.
LoglikWeibullSurvRegCens(x, data_y, data_delta_loglik, data_cov_noncens = NULL, data_cov_cens, density, data_r_loglik, data_lowerbound, intlimit = 10^-10)
LoglikWeibullSurvRegCens(x, data_y, data_delta_loglik, data_cov_noncens = NULL, data_cov_cens, density, data_r_loglik, data_lowerbound, intlimit = 10^-10)
x |
Vector of parameters, ordered as follows: Scale parameter, Shape parameter, regression parameters (i.e. |
data_y |
Time-to-event vector. |
data_delta_loglik |
Censored indicator vector of the time-to-event (0: censored, 1: not censored). |
data_cov_noncens |
Matrix where each column represents a non-censored covariate. |
data_cov_cens |
Censored covariate vector. |
density |
Density function of the censored covariate. |
data_r_loglik |
Censored indicator vector of the censored covariate (0: censored, 1: not censored). |
data_lowerbound |
A vector which corresponds to the lower bounds for the interval-censored observations of the censored covariate. If no lower bound is available then put |
intlimit |
In computation of integrals, values of the function to be integrated below |
Function not intended to be invoked by the user.
Stanislas Hubeaux, [email protected]
Kaspar Rufibach, [email protected]
http://www.kasparrufibach.ch
Hubeaux, S. (2013). Parametric Surival Regression Model with left- and/or interval-censored covariate. Technical report, Biostatistics Oncology, F. Hoffmann-La Roche Ltd.
Sattar, A., Sinha, S. K. and Morris, N. J. (2012). A Parametric Survival Model When a Covariate is Subject to Left-Censoring. Biometrics & Biostatistics, S3(2).
Computes estimates of the parameters of two censored Normal samples, as well as the mean difference between the two samples.
NormalMeanDiffCens(censdata1, censdata2, conf.level = 0.95, null.values = c(0, 0, 1, 1))
NormalMeanDiffCens(censdata1, censdata2, conf.level = 0.95, null.values = c(0, 0, 1, 1))
censdata1 |
Observations of first sample, format as specified by |
censdata2 |
Observations of second sample, as specified by |
conf.level |
Confidence level for confidence intervals. |
null.values |
Fixed values for hypothesis tests. Four-dimensional vector specifying the hypothesis for |
A table with estimators and inference for the means and standard deviations of both samples, as well as the difference between the mean of the first and second sample. Hypothesis tests are for the values in
null.values
and for the null hypothesis of no mean difference.
Stanislas Hubeaux, [email protected]
Kaspar Rufibach, [email protected]
http://www.kasparrufibach.ch
Hubeaux, S. (2013). Estimation from left- and/or interval-censored samples. Technical report, Biostatistics Oncology, F. Hoffmann-La Roche Ltd.
Lynn, H. S. (2001). Maximum likelihood inference for left-censored HIV RNA data. Stat. Med., 20, 33–45.
## example with interval-censored Normal samples n <- 500 prop.cens <- 0.35 mu <- c(0, 2) sigma <- c(1, 1) set.seed(2013) ## Sample 1: LOD1 <- qnorm(prop.cens, mean = mu[1], sd = sigma[1]) x1 <- rnorm(n, mean = mu[1], sd = sigma[1]) s1 <- censorContVar(x1, LLOD = LOD1) ## Sample 2: LOD2 <- qnorm(0.35, mean = mu[2], sd = sigma[2]) x2 <- rnorm(n, mean = mu[2], sd = sigma[2]) s2 <- censorContVar(x2, LLOD = LOD2) ## inference on distribution parameters and mean difference: NormalMeanDiffCens(censdata1 = s1, censdata2 = s2)
## example with interval-censored Normal samples n <- 500 prop.cens <- 0.35 mu <- c(0, 2) sigma <- c(1, 1) set.seed(2013) ## Sample 1: LOD1 <- qnorm(prop.cens, mean = mu[1], sd = sigma[1]) x1 <- rnorm(n, mean = mu[1], sd = sigma[1]) s1 <- censorContVar(x1, LLOD = LOD1) ## Sample 2: LOD2 <- qnorm(0.35, mean = mu[2], sd = sigma[2]) x2 <- rnorm(n, mean = mu[2], sd = sigma[2]) s2 <- censorContVar(x2, LLOD = LOD2) ## inference on distribution parameters and mean difference: NormalMeanDiffCens(censdata1 = s1, censdata2 = s2)
Computes maximum likelihood estimators of the canonical parameters for several distributions, based on a censored sample.
ParamSampleCens(censdata, dist = c("normal", "logistic", "gamma", "weibull")[1], null.values = c(0, 1), conf.level = 0.95, initial = NULL)
ParamSampleCens(censdata, dist = c("normal", "logistic", "gamma", "weibull")[1], null.values = c(0, 1), conf.level = 0.95, initial = NULL)
censdata |
Dataframe that contains censored data, format as specified by |
dist |
Assumed distribution of the sample. |
null.values |
Fixed values for hypothesis tests. |
conf.level |
Confidence level of confidence intervals. |
initial |
Initial values for the maximization. |
coeff |
Estimators, standard errors, confidence intervals, and 2-sided |
percent.cens |
Percentage of censored observations. |
loglik |
Log likelihood function value at the estimator. |
info.converg |
Convergence information provided by the function |
info.converg.message |
Message provided by the function |
Functions with similar functionality are provided in the package fitdistrplus.
Stanislas Hubeaux, [email protected]
Kaspar Rufibach, [email protected]
http://www.kasparrufibach.ch
Hubeaux, S. (2013). Estimation from left- and/or interval-censored samples. Technical report, Biostatistics Oncology, F. Hoffmann-La Roche Ltd.
Lynn, H. S. (2001). Maximum likelihood inference for left-censored HIV RNA data. Stat. Med., 20, 33–45.
n <- 500 prop.cens <- 0.35 ## example with a left-censored Normally distributed sample set.seed(2013) mu <- 3.5 sigma <- 1 LOD <- qnorm(prop.cens, mean = mu, sd = sigma) x1 <- rnorm(n, mean = mu, sd = sigma) s1 <- censorContVar(x1, LLOD = LOD) ParamSampleCens(censdata = s1) ## example with an interval-censored Normal sample set.seed(2013) x2 <- rnorm(n, mean = mu, sd = sigma) LOD <- qnorm(prop.cens / 2, mean = mu, sd = sigma) UOD <- qnorm(1 - prop.cens / 2, mean = mu, sd = sigma) s2 <- censorContVar(x2, LLOD = LOD, ULOD = UOD) ParamSampleCens(censdata = s2) ## Not run: ## compare to fitdistrplus library(fitdistrplus) s2 <- as.data.frame(s2) colnames(s2) <- c("left", "right") summary(fitdistcens(censdata = s2, distr = "norm")) ## End(Not run)
n <- 500 prop.cens <- 0.35 ## example with a left-censored Normally distributed sample set.seed(2013) mu <- 3.5 sigma <- 1 LOD <- qnorm(prop.cens, mean = mu, sd = sigma) x1 <- rnorm(n, mean = mu, sd = sigma) s1 <- censorContVar(x1, LLOD = LOD) ParamSampleCens(censdata = s1) ## example with an interval-censored Normal sample set.seed(2013) x2 <- rnorm(n, mean = mu, sd = sigma) LOD <- qnorm(prop.cens / 2, mean = mu, sd = sigma) UOD <- qnorm(1 - prop.cens / 2, mean = mu, sd = sigma) s2 <- censorContVar(x2, LLOD = LOD, ULOD = UOD) ParamSampleCens(censdata = s2) ## Not run: ## compare to fitdistrplus library(fitdistrplus) s2 <- as.data.frame(s2) colnames(s2) <- c("left", "right") summary(fitdistcens(censdata = s2, distr = "norm")) ## End(Not run)
print
method for class "src"
.
## S3 method for class 'src' print(x, ...)
## S3 method for class 'src' print(x, ...)
x |
An object of class |
... |
Further arguments. |
The function print.src
returns the estimated parameters of the Weibull regression, incl. AIC, when calling SurvRegCens
.
Stanislas Hubeaux, [email protected]
Kaspar Rufibach, [email protected]
http://www.kasparrufibach.ch
Hubeaux, S. (2013). Parametric Surival Regression Model with left- and/or interval-censored covariate. Technical report, Biostatistics Oncology, F. Hoffmann-La Roche Ltd.
Hubeaux, S. and Rufibach, K. (2014). SurvRegCensCov: Weibull Regression for a Right-Censored Endpoint with a Censored Covariate. Preprint, https://arxiv.org/abs/1402.0432.
Sattar, A., Sinha, S. K. and Morris, N. J. (2012). A Parametric Survival Model When a Covariate is Subject to Left-Censoring. Biometrics & Biostatistics, S3(2).
## See help file of function "SurvRegCens".
## See help file of function "SurvRegCens".
summary
method for class "src"
.
## S3 method for class 'src' summary(object, ...)
## S3 method for class 'src' summary(object, ...)
object |
An object of class |
... |
Further arguments. |
The function summary.src
returns the estimated parameters, incl. statistical inference, of the Weibull regression, incl. AIC, when calling SurvRegCens
.
Stanislas Hubeaux, [email protected]
Kaspar Rufibach, [email protected]
http://www.kasparrufibach.ch
Hubeaux, S. (2013). Parametric Surival Regression Model with left- and/or interval-censored covariate. Technical report, Biostatistics Oncology, F. Hoffmann-La Roche Ltd.
Hubeaux, S. and Rufibach, K. (2014). SurvRegCensCov: Weibull Regression for a Right-Censored Endpoint with a Censored Covariate. Preprint, https://arxiv.org/abs/1402.0432.
Sattar, A., Sinha, S. K. and Morris, N. J. (2012). A Parametric Survival Model When a Covariate is Subject to Left-Censoring. Biometrics & Biostatistics, S3(2).
## See help file of function "SurvRegCens".
## See help file of function "SurvRegCens".
Computes estimators for the shape and scale parameter of the Weibull distribution, as well as for the vector of regression parameters in a parametric survival model with potentially right-censored time-to-event endpoint distributed according to a Weibull distribution. The regression allows for one potentially interval-censored and an arbitrary number of non-censored covariates.
SurvRegCens(formula, data = parent.frame(), Density, initial, conf.level = 0.95, intlimit = 10^-10, namCens = "VarCens", trace = 0, reltol = 10^-8)
SurvRegCens(formula, data = parent.frame(), Density, initial, conf.level = 0.95, intlimit = 10^-10, namCens = "VarCens", trace = 0, reltol = 10^-8)
formula |
A formula expression as for other regression models. The response has to be a survival object for right-censored data, as returned by the |
data |
A data frame in which to interpret the variables named in the formula argument. |
Density |
Density function of the censored covariate. |
initial |
Initial values for the parameters to be optimized over, ordered according to Scale parameter, Shape parameter, regression parameters (i.e. |
conf.level |
Confidence level of confidence intervals. |
intlimit |
In computation of integrals, values of the function to be integrated below |
namCens |
Name of censored covariate, to tidy outputs. |
trace |
|
reltol |
|
The time-to-event distributed according to a Weibull distribution, i.e. time-to-event Weibull
, has conditional density given by,
conditional hazard function given by,
and conditional survival function given by,
where collects the values of each covariate for observation
and
represents the regression parameters.
SurvRegCens
returns an object of class "src"
, a list containing the
following components:
coeff |
Estimators, confidence intervals, |
percent.cens |
Percentage of censored observations in the censored covariate. |
loglik |
Log-likelihood function value at the estimators. |
info.converg |
Convergence information provided by the function |
info.converg.message |
Message provided by |
The methods print.src
, summary.src
, coef.src
, and logLik.src
are used to print or obtain a summary, coefficients, or the value of the log-likelihood at the maximum.
Stanislas Hubeaux, [email protected]
Kaspar Rufibach, [email protected]
http://www.kasparrufibach.ch
Hubeaux, S. (2013). Parametric Surival Regression Model with left- and/or interval-censored covariate. Technical report, Biostatistics Oncology, F. Hoffmann-La Roche Ltd.
Hubeaux, S. and Rufibach, K. (2014). SurvRegCensCov: Weibull Regression for a Right-Censored Endpoint with a Censored Covariate. Preprint, https://arxiv.org/abs/1402.0432.
Sattar, A., Sinha, S. K. and Morris, N. J. (2012). A Parametric Survival Model When a Covariate is Subject to Left-Censoring. Biometrics & Biostatistics, S3(2).
## Not run: ## -------------------------------------------------------------- ## 1 censored-covariate and 2 non-censored covariates ## no censoring, to compare result with survival::survreg ## modify prop.cens to introduce left-censoring of covariate ## -------------------------------------------------------------- set.seed(158) n <- 100 lambda <- exp(-2) gamma <- 1.5 ## vector of regression parameters: the last entry is the one for the censored covariate beta <- c(0.3, -0.2, 0.25) true <- c(lambda, gamma, beta) ## non-censored covariates var1 <- rnorm(n, mean = 4, sd = 0.5) var2 <- rnorm(n, mean = 4, sd = 0.5) ## Generate censored covariate. ## For generation of Weibull survival times, do not left-censor it yet. var3 <- rnorm(n, mean = 5, sd = 0.5) ## simulate from a Weibull regression model time <- TimeSampleWeibull(covariate_noncens = data.frame(var1, var2), covariate_cens = var3, lambda = lambda, gamma = gamma, beta = beta) ## left-censor covariate ## prop.cens specifies the proportion of observations that should be left-censored prop.cens <- 0 LOD <- qnorm(prop.cens, mean = 5, sd = 0.5) var3.cens <- censorContVar(var3, LLOD = LOD) ## censor survival time event <- matrix(1, nrow = n, ncol = 1) time.cens <- rexp(n, rate = 0.5) ind.time <- (event >= time.cens) event[ind.time] <- 0 time[ind.time] <- time.cens[ind.time] ## specify the density for the censored covariate: ## For simplicity, we take here the "true" density we simulate from. In an application, ## you might want to use a density with parameters estimated from the censored covariate, ## e.g. using the function ParamSampleCens. See example in Hubeaux & Rufibach (2014). DensityCens <- function(value){return(dnorm(value, mean = 5, sd = 0.5))} ## use Weibull regression where each censored covariate value is set ## to LOD ("naive" method) naive <- survreg(Surv(time, event) ~ var1 + var2 + var3.cens[, 2], dist = "weibull") initial <- as.vector(ConvertWeibull(naive)$vars[, 1]) ## use new method that takes into account the left-censoring of one covariate data <- data.frame(time, event, var3.cens, var1, var2) formula <- formula(Surv(time, event) ~ Surv(time = var3.cens[, 1], time2 = var3.cens[, 2], type = "interval2") + var1 + var2) cens1 <- SurvRegCens(formula = formula, data = data, Density = DensityCens, initial = initial, namCens = "biomarker") summary(cens1) coef(cens1) logLik(cens1) ## compare estimates tab <- data.frame(cbind(true, initial, cens1$coeff[, 1])) colnames(tab) <- c("true", "naive", "Weibull MLE") rownames(tab) <- rownames(cens1$coeff) tab ## compare confidence intervals ConvertWeibull(naive)$HR[, 2:3] cens1$coeff[, 7:8] ## -------------------------------------------------------------- ## model without the non-censored covariates ## -------------------------------------------------------------- naive2 <- survreg(Surv(time, event) ~ var3.cens[, 2], dist = "weibull") initial2 <- as.vector(ConvertWeibull(naive2)$vars[, 1]) ## use new method that takes into account the left-censoring of one covariate formula <- formula(Surv(time, event) ~ Surv(time = var3.cens[, 1], time2 = var3.cens[, 2], type = "interval2")) cens2 <- SurvRegCens(formula = formula, data = data, Density = DensityCens, initial = initial2, namCens = "biomarker") summary(cens2) ## compare estimates tab <- data.frame(cbind(true[c(1, 2, 5)], initial2, cens2$coeff[, 1])) colnames(tab) <- c("true", "naive", "Weibull MLE") rownames(tab) <- rownames(cens2$coeff) tab ## compare confidence intervals ConvertWeibull(naive2)$HR[, 2:3] cens2$coeff[, 7:8] ## End(Not run)
## Not run: ## -------------------------------------------------------------- ## 1 censored-covariate and 2 non-censored covariates ## no censoring, to compare result with survival::survreg ## modify prop.cens to introduce left-censoring of covariate ## -------------------------------------------------------------- set.seed(158) n <- 100 lambda <- exp(-2) gamma <- 1.5 ## vector of regression parameters: the last entry is the one for the censored covariate beta <- c(0.3, -0.2, 0.25) true <- c(lambda, gamma, beta) ## non-censored covariates var1 <- rnorm(n, mean = 4, sd = 0.5) var2 <- rnorm(n, mean = 4, sd = 0.5) ## Generate censored covariate. ## For generation of Weibull survival times, do not left-censor it yet. var3 <- rnorm(n, mean = 5, sd = 0.5) ## simulate from a Weibull regression model time <- TimeSampleWeibull(covariate_noncens = data.frame(var1, var2), covariate_cens = var3, lambda = lambda, gamma = gamma, beta = beta) ## left-censor covariate ## prop.cens specifies the proportion of observations that should be left-censored prop.cens <- 0 LOD <- qnorm(prop.cens, mean = 5, sd = 0.5) var3.cens <- censorContVar(var3, LLOD = LOD) ## censor survival time event <- matrix(1, nrow = n, ncol = 1) time.cens <- rexp(n, rate = 0.5) ind.time <- (event >= time.cens) event[ind.time] <- 0 time[ind.time] <- time.cens[ind.time] ## specify the density for the censored covariate: ## For simplicity, we take here the "true" density we simulate from. In an application, ## you might want to use a density with parameters estimated from the censored covariate, ## e.g. using the function ParamSampleCens. See example in Hubeaux & Rufibach (2014). DensityCens <- function(value){return(dnorm(value, mean = 5, sd = 0.5))} ## use Weibull regression where each censored covariate value is set ## to LOD ("naive" method) naive <- survreg(Surv(time, event) ~ var1 + var2 + var3.cens[, 2], dist = "weibull") initial <- as.vector(ConvertWeibull(naive)$vars[, 1]) ## use new method that takes into account the left-censoring of one covariate data <- data.frame(time, event, var3.cens, var1, var2) formula <- formula(Surv(time, event) ~ Surv(time = var3.cens[, 1], time2 = var3.cens[, 2], type = "interval2") + var1 + var2) cens1 <- SurvRegCens(formula = formula, data = data, Density = DensityCens, initial = initial, namCens = "biomarker") summary(cens1) coef(cens1) logLik(cens1) ## compare estimates tab <- data.frame(cbind(true, initial, cens1$coeff[, 1])) colnames(tab) <- c("true", "naive", "Weibull MLE") rownames(tab) <- rownames(cens1$coeff) tab ## compare confidence intervals ConvertWeibull(naive)$HR[, 2:3] cens1$coeff[, 7:8] ## -------------------------------------------------------------- ## model without the non-censored covariates ## -------------------------------------------------------------- naive2 <- survreg(Surv(time, event) ~ var3.cens[, 2], dist = "weibull") initial2 <- as.vector(ConvertWeibull(naive2)$vars[, 1]) ## use new method that takes into account the left-censoring of one covariate formula <- formula(Surv(time, event) ~ Surv(time = var3.cens[, 1], time2 = var3.cens[, 2], type = "interval2")) cens2 <- SurvRegCens(formula = formula, data = data, Density = DensityCens, initial = initial2, namCens = "biomarker") summary(cens2) ## compare estimates tab <- data.frame(cbind(true[c(1, 2, 5)], initial2, cens2$coeff[, 1])) colnames(tab) <- c("true", "naive", "Weibull MLE") rownames(tab) <- rownames(cens2$coeff) tab ## compare confidence intervals ConvertWeibull(naive2)$HR[, 2:3] cens2$coeff[, 7:8] ## End(Not run)
Generates time-to-event data using the transform inverse sampling method, and such that the time-to-event is distributed according to a Weibull distribution induced by censored and/or non-censored covariates. Can be used to set up simulations.
TimeSampleWeibull(covariate_noncens = NULL, covariate_cens, lambda, gamma, beta)
TimeSampleWeibull(covariate_noncens = NULL, covariate_cens, lambda, gamma, beta)
covariate_cens |
Censored covariate vector. |
covariate_noncens |
Matrix where each column represents a non-censored covariate. |
lambda |
Scale parameter. |
gamma |
Shape parameter. |
beta |
Regression parameters, ordered as |
The use of this function is illustrated in SurvRegCens
.
Stanislas Hubeaux, [email protected]
This function constructs a diagnostic plot of the adequacy of the Weibull distribution for survival data with respect to one categorical covariate. If the Weibull distribution fits the data well, then the lines produced should be linear and parallel.
WeibullDiag(formula, data = parent.frame(), labels = names(m$strata))
WeibullDiag(formula, data = parent.frame(), labels = names(m$strata))
formula |
A formula containing a |
data |
Data set. |
labels |
A vector containing labels for the plotted lines. |
As discussed in Klein and Moeschberger (2003), one method for checking the adequacy of the Weibull model with a
categorical covariate is to produce stratified Kaplan-Meier estimates (KM), which can be transformed to estimate
the log cumulative hazard for each stratum. Then in a plot of versus
, the
lines should be linear and parallel. This can be seen as the log cumulative hazard for the Weibull distribution
is
Produces a plot of log Time vs. log Estimated Cumulative Hazard for each level of the predictor
(similarly to what can be obtained using plot.survfit
and the fun = "cloglog"
option),
as well as a data set containing that information.
Sarah R. Haile, Epidemiology, Biostatistics and Prevention Institute (EBPI), University of Zurich, [email protected]
Klein, J. and Moeschberger, M. (2003). Survival analysis: techniques for censored and truncated data. 2nd edition, Springer.
Requires survival. A similar plot can be produced using plot.survfit
and the option fun = "cloglog"
.
data(larynx) WeibullDiag(Surv(time, death) ~ stage, data = larynx)
data(larynx) WeibullDiag(Surv(time, death) ~ stage, data = larynx)
SurvRegCens
Function to be integrated to compute log-likelihood function for the Weibull survival regression model with a censored covariate.
WeibullIntegrate(x, x_i_noncens = NULL, density, param_y_i, param_delta_i, param_lambda, param_gamma, param_beta, intlimit = 10^-10, ForIntegrate = TRUE)
WeibullIntegrate(x, x_i_noncens = NULL, density, param_y_i, param_delta_i, param_lambda, param_gamma, param_beta, intlimit = 10^-10, ForIntegrate = TRUE)
x |
Value of the censored covariate for observation |
x_i_noncens |
Vector of values of the non-censored covariates for observation |
density |
Density function of the censored covariate. |
param_y_i |
Value of the time-to-event for observation |
param_delta_i |
Censoring indicator of time-to-event for observation |
param_lambda |
Scale parameter of the Weibull distribution. |
param_gamma |
Shape parameter of the Weibull distribution. |
param_beta |
Regression parameters (i.e. |
intlimit |
In computation of integrals, values of the function to be integrated below |
ForIntegrate |
|
Function is not intended to be invoked by the user.
Stanislas Hubeaux, [email protected]
WeibullReg
performs Weibull regression using the survreg
function, and transforms the
estimates to a more natural parameterization. Additionally, it produces hazard ratios (corresponding to the proportional
hazards interpretation), and event time ratios (corresponding to the accelerated failure time interpretation) for all
covariates.
WeibullReg(formula, data = parent.frame(), conf.level = 0.95)
WeibullReg(formula, data = parent.frame(), conf.level = 0.95)
formula |
A |
data |
The dataset containing all variables referenced in |
conf.level |
Specifies that |
Details regarding the transformations of the parameters and their standard errors can be found in Klein and
Moeschberger (2003, Chapter 12). An explanation of event time ratios for the accelerated failure time interpretation of
the model can be found in Carroll (2003). A general overview can be found in the vignette("weibull")
of this package, or in the documentation for ConvertWeibull
.
formula |
The formula for the Weibull regression model. |
coef |
The transformed maximum likelihood estimates, with standard errors. |
HR |
The hazard ratios for each of the predictors, with |
ETR |
The event time ratios (acceleration factors) for each of the predictors, with |
summary |
The summary output from the original |
Sarah R. Haile, Epidemiology, Biostatistics and Prevention Institute (EBPI), University of Zurich, [email protected]
Carroll, K. (2003). On the use and utility of the Weibull model in the analysis of survival data. Controlled Clinical Trials, 24, 682–701.
Klein, J. and Moeschberger, M. (2003). Survival analysis: techniques for censored and truncated data. 2nd edition, Springer.
Requires the package survival. This function depends on ConvertWeibull
. See also
survreg
.
data(larynx) WR <- WeibullReg(Surv(time, death) ~ factor(stage) + age, data = larynx) WR
data(larynx) WR <- WeibullReg(Surv(time, death) ~ factor(stage) + age, data = larynx) WR