Description Usage Arguments Details Value Examples
View source: R/impute2CoxSurv.R
Performs survival analysis using Cox proportional hazard models on imputed genetic data from IMPUTE2 output
1 2 3 4 5 6 | impute2CoxSurv(impute.file, sample.file, chr, covariate.file, id.column,
sample.ids = NULL, time.to.event, event, covariates,
inter.term = NULL, print.covs = "only", out.file,
chunk.size = 10000, maf.filter = 0.05, exclude.snps = NULL,
flip.dosage = TRUE, verbose = TRUE, clusterObj = NULL,
keepGDS = FALSE)
|
impute.file |
character of IMPUTE2 file |
sample.file |
character of sample file affiliated with IMPUTE2 file |
chr |
numeric denoting chromosome number |
covariate.file |
data.frame comprising phenotype information, all covariates to be added in the model must be numeric. |
id.column |
character giving the name of the ID column in covariate.file. |
sample.ids |
character vector of sample IDs to keep in survival analysis |
time.to.event |
character of column name in covariate.file that represents the time interval of interest in the analysis |
event |
character of column name in covariate.file that represents the event of interest to be included in the analysis |
covariates |
character vector with exact names of columns in covariate.file to include in analysis |
inter.term |
character string giving the column name of the covariate
that will be added to the interaction term with
SNP (e.g. |
print.covs |
character string of either |
out.file |
character of output file name (do not include extension) |
chunk.size |
integer number of variants to process per thread |
maf.filter |
numeric to filter minor allele frequency
(i.e. choosing 0.05 means filtering MAF>0.05). User can set this to
|
exclude.snps |
a character vector listing the rsIDs of SNPs that will be excluded from analyses |
flip.dosage |
logical to flip which allele the dosage was
calculated on, default |
verbose |
logical for messages that describe which part of the analysis is currently being run |
clusterObj |
A cluster object that can be used with the
|
keepGDS |
logical to keep GDS files (compressed IMPUTE2 files)
after the analysis. Defaults to |
Testing for SNP-covariate interactions:
User can define the column name of the covariate that will be included in
the interaction term.
For example, for given covariates a
and b
, where c
is
defined as the inter.term
the model will be:
~ a + b + c + SNP + c*SNP
.
Printing results of other covariates:
print.covs
argument controls the number of covariates will be printed
as output. The function is set to only
by default and will only print the SNP or if an interaction term is given,
the results of the interaction
term (e.g. SNP*covariate
). Whereas, all
will print results
(coef, se.coef, p.value etc) of all covariates included in the model.
some
is only applicable if an interaction term is given and will
print the results for SNP, covariate tested for interaction and the
interaction term. User should be mindful about using the all
option,
as it will likely slow down the analysis and will increase
the output file size.
User defined parallelization:
This function uses parApply
from parallel
package to fit
models to SNPs in parallel.
User is not required to set any options for the parallelization.
However, advanced users who wish to optimize it, can provide a
cluster object generated by makeCluster
family of functions that
suits their need and platform.
Saves two text files directly to disk:
.coxph
extension containing CoxPH survival analysis results.
.snps_removed
extension containing SNPs that were removed
due to low variance or user-defined thresholds.
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 | impute.file <- system.file(package="gwasurvivr",
"extdata",
"impute_example.impute2.gz")
sample.file <- system.file(package="gwasurvivr",
"extdata",
"impute_example.impute2_sample")
covariate.file <- system.file(package="gwasurvivr",
"extdata",
"simulated_pheno.txt")
covariate.file <- read.table(covariate.file,
sep=" ",
header=TRUE,
stringsAsFactors = FALSE)
covariate.file$SexFemale <- ifelse(covariate.file$sex=="female", 1L, 0L)
sample.ids <- covariate.file[covariate.file$group=="experimental",]$ID_2
impute2CoxSurv(impute.file=impute.file,
sample.file=sample.file,
chr=14,
covariate.file=covariate.file,
id.column="ID_2",
sample.ids=sample.ids,
time.to.event="time",
event="event",
covariates=c("age", "SexFemale", "DrugTxYes"),
inter.term=NULL,
print.covs="only",
out.file="impute_example",
chunk.size=50,
maf.filter=0.005,
exclude.snps=NULL,
flip.dosage=TRUE,
verbose=TRUE,
clusterObj=NULL,
keepGDS=FALSE)
|
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.