gdsCoxSurv: Fit Cox survival to all variants from a standard IMPUTE2...

Description Usage Arguments Details Value Examples

View source: R/gdsCoxSurv.R

Description

Performs survival analysis using Cox proportional hazard models on imputed genetic data from GDS files.

Usage

1
2
3
4
5
gdsCoxSurv(gdsfile, covariate.file, id.column, sample.ids = NULL,
  time.to.event, event, covariates, inter.term = NULL,
  print.covs = "only", out.file, chunk.size = 5000,
  maf.filter = 0.05, flip.dosage = TRUE, verbose = TRUE,
  clusterObj = NULL)

Arguments

gdsfile

path to .gds file. Location of the .gds file should also contain .snp.rdata and .scan.rdata files.

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. term*SNP). See details.

print.covs

character string of either "only", "all" or "some", defining which covariate statistics should be printed to the output. See details.

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 NULL if no filtering is preffered. Default is 0.05.

flip.dosage

logical to flip which allele the dosage was calculated on, default flip.dosage=TRUE

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 parApply function. See details.

Details

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.

Value

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.

Examples

 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
gdsfile <- system.file(package="gwasurvivr",
                       "extdata",
                       "gds_example.gds")
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
gdsCoxSurv(gdsfile=gdsfile,
           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,
           flip.dosage=TRUE,
           verbose=TRUE,
           clusterObj=NULL)  
 

gwasurvivr documentation built on Nov. 8, 2020, 6:53 p.m.