.
.
options(width=110)
A phenome-wide association study (PheWAS) is known to be a powerful tool in discovery and replication of genetic association studies. The recent development in the UK Biobank resource with deep genomic and phenotyping data has provided unparalleled research opportunities. To reduce the computational complexity and cost of PheWAS, the SAIGE (scalable and accurate implementation of generalized mixed model [1]) method was proposed recently, controlling for case-control imbalance and sample structure in single variant association studies. However, it is still computationally challenging to analyze the associations of thousands of phenotypes with whole-genome variant data, especially for disease diagnoses using the ICD-10 codes of deep phenotyping.
Here we develop a new high-performance statistical package (SAIGEgds) for large-scale PheWAS using mixed models [2]. In this package, the SAIGE method is implemented with optimized C++ codes taking advantage of sparse structure of genotype dosages. SAIGEgds supports efficient genomic data structure (GDS) files [3] including both integer genotypes and numeric imputed dosages. Benchmarks using the UK Biobank White British genotype data (N=430K) with coronary heart disease and simulated cases, show that SAIGEgds is 5 to 6 times faster than the SAIGE R package in the steps of fitting null models and p-value calculations. When used in conjunction with high-performance computing (HPC) clusters and/or cloud resources, SAIGEgds provides an efficient analysis pipeline for biobank-scale PheWAS.
.
if (!requireNamespace("BiocManager", quietly=TRUE)) install.packages("BiocManager") BiocManager::install("SAIGEgds")
The BiocManager::install()
approach may require that you build from source, i.e. make
and compilers must be installed on your system -- see the R FAQ for your operating system; you may also need to install dependencies manually.
.
.
library(SeqArray) library(SAIGEgds) # the genotype file in SeqArray GDS format (1000 Genomes Phase1, chromosome 22) (geno_fn <- seqExampleFileName("KG_Phase1"))
# open a SeqArray file in the package (1000 Genomes Phase1, chromosome 22) gds <- seqOpen(geno_fn)
The LD pruning is provided by snpgdsLDpruning() in the pacakge SNPRelate:
library(SNPRelate) set.seed(1000) snpset <- snpgdsLDpruning(gds) str(snpset) snpset.id <- unlist(snpset, use.names=FALSE) # get the variant IDs of a LD-pruned set head(snpset.id)
Create a genotype file for genetic relationship matrix (GRM) using the LD-pruned SNP set:
grm_fn <- "grm_geno.gds" seqSetFilter(gds, variant.id=snpset.id) # export to a GDS genotype file without annotation data seqExport(gds, grm_fn, info.var=character(), fmt.var=character(), samp.var=character())
If genotypes are split by chromosomes, seqMerge()
in the SeqArray package can be used to combine the LD-pruned SNP sets.
# close the file seqClose(gds)
.
A simulated phenotype data is used to demonstrate the model fitting:
set.seed(1000) sampid <- seqGetData(grm_fn, "sample.id") # sample IDs in the genotype file pheno <- data.frame( sample.id = sampid, y = sample(c(0, 1), length(sampid), replace=TRUE, prob=c(0.95, 0.05)), x1 = rnorm(length(sampid)), x2 = rnorm(length(sampid)), stringsAsFactors = FALSE) head(pheno) # null model fitting using GRM from grm_fn grm_fn glmm <- seqFitNullGLMM_SPA(y ~ x1 + x2, pheno, grm_fn, trait.type="binary", sample.col="sample.id")
.
# genetic variants stored in the file geno_fn geno_fn # calculate, using 2 processes assoc <- seqAssocGLMM_SPA(geno_fn, glmm, mac=10, parallel=2) head(assoc) # filtering based on p-value assoc[assoc$pval < 5e-4, ]
The output could be directly saved to an R object file or a GDS file:
# save to 'assoc.gds' seqAssocGLMM_SPA(geno_fn, glmm, mac=10, parallel=2, res.savefn="assoc.gds")
Open the output GDS file using the functions in the gdsfmt package:
# open the GDS file (f <- openfn.gds("assoc.gds")) # get p-values pval <- read.gdsn(index.gdsn(f, "pval")) summary(pval) closefn.gds(f)
Load association results using the function seqSAIGE_LoadPval()
in SAIGEgds:
res <- seqSAIGE_LoadPval("assoc.gds") head(res)
.
n <- nrow(assoc) # expected p-values exp.pval <- (1:n) / n # observed p-values obs.pval <- sort(assoc$pval) plot(-log10(exp.pval), -log10(obs.pval), xlab="-log10(expected P)", ylab="-log10(observed P)", cex=2/3) abline(0, 1, col="red", lty=2)
.
.
sessionInfo()
unlink(c("grm_geno.gds", "assoc.gds"), force=TRUE)
.
.
SeqArray: Data Management of Large-scale Whole-genome Sequence Variant Calls
SNPRelate: Parallel Computing Toolset for Relatedness and Principal Component Analysis of SNP Data
gds2bgen: Format conversion from bgen to gds
Any scripts or data that you put into this service are public.
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.