assocTestAggregate | R Documentation |
assocTestAggregate
performs aggregate association tests using the null model fit with fitNullModel
.
## S4 method for signature 'SeqVarIterator'
assocTestAggregate(gdsobj, null.model, AF.max=1,
weight.beta=c(1,1), weight.user=NULL,
test=c("Burden", "SKAT", "fastSKAT", "SMMAT", "fastSMMAT",
"SKATO", "BinomiRare", "CMP"),
neig = 200, ntrace = 500,
rho = seq(from = 0, to = 1, by = 0.1),
sparse=TRUE, imputed=FALSE,
male.diploid=TRUE, genome.build=c("hg19", "hg38"),
BPPARAM=bpparam(), verbose=TRUE)
## S4 method for signature 'GenotypeIterator'
assocTestAggregate(gdsobj, null.model, AF.max=1,
weight.beta=c(1,1), weight.user=NULL,
test=c("Burden", "SKAT", "fastSKAT", "SMMAT", "fastSMMAT",
"SKATO", "BinomiRare", "CMP"),
neig = 200, ntrace = 500,
rho = seq(from = 0, to = 1, by = 0.1),
male.diploid=TRUE, BPPARAM=bpparam(), verbose=TRUE)
gdsobj |
An object of class |
null.model |
A null model object returned by |
AF.max |
A numeric value specifying the upper bound on the effect allele frequency for variants to be included in the analysis. |
weight.beta |
A numeric vector of length two specifying the two parameters of the Beta distribution used to determine variant weights; weights are given by |
weight.user |
A character string specifying the name of a variable to be used as variant weights. This variable can be in either 1) the variantData slot of |
test |
A character string specifying the type of test to be performed. The possibilities are |
neig |
The number eigenvalues to approximate by using random projections for calculating p-values with fastSKAT or fastSMMAT; default is 200. See 'Details' for more information. |
ntrace |
The number of vectors to sample when using random projections to estimate the trace needed for p-value calculation with fastSKAT or fastSMMAT; default is 500. See 'Details' for more information. |
rho |
A numeric value (or vector of numeric values) in |
sparse |
Logical indicator of whether to read genotypes as sparse Matrix objects; the default is |
imputed |
Logical indicator of whether to read dosages from the "DS" field containing imputed dosages instead of counting the number of alternate alleles. |
male.diploid |
Logical for whether males on sex chromosomes are coded as diploid. Default is 'male.diploid=TRUE', meaning sex chromosome genotypes for males have values 0/2. If the input |
genome.build |
A character sting indicating genome build; used to identify pseudoautosomal regions on the X and Y chromosomes. These regions are not treated as sex chromosomes when calculating allele frequencies. |
BPPARAM |
A |
verbose |
Logical indicator of whether updates from the function should be printed to the console; the default is |
The type of aggregate unit tested depends on the class of iterator
used for gdsobj
. Options include sliding windows, specific
ranges of variants or selection of individual variants (ranges with
width 1). See SeqVarIterator
for more details.
assocTestAggregate
uses the BiocParallel
package to process iterator chunks in parallel. See the BiocParallel
documentation for more information on the default behaviour of bpparam
and how to register different parallel backends. If serial execution is desired, set BPPARAM=BiocParallel::SerialParam()
. Note that parallel execution requires more RAM than serial execution.
All samples included in null model
will be included in the
association test, even if a different set of samples is present in the
current filter for gdsobj
.
The effect size estimate is for each copy of the alternate allele (when gdsobj
is a SeqVarIterator
object) or the "A" allele (when gdsobj
is a GenotypeIterator
object). We refer to this as the "effect allele" in the rest of the documentation. For multiallelic variants in SeqVarIterator
objects, each alternate (or "A") allele is tested separately.
Monomorphic variants (including variants where every sample is a heterozygote) are always omitted from the aggregate unit prior to testing.
Somewhat similarly to SKAT-O, the variant Set Mixed Model Association Test (SMMAT, Chen et al., 2019) combines the burden test p-value with an adjusted SKAT (which is asymptotically independent of the burden test) p-value using a chi-square distribution with 4df from Fisher's method.
SKAT and SMMAT will attempt to use Davies' method (i.e. integration) to calculate p-values; if an error occurs in integration or the reported p-values are too small that they are unreliable (i.e. near machine epsilon), then the saddlepoint approximation will instead be used to calculate the p-values.
The fastSKAT method of Lumley et al. (2018) uses random matrix theory to speed up the computation of SKAT p-values. When test = "fastSKAT"
, the function attempts to inteligently determine which p-value calculation approach to use for each aggregation unit: (1) if min(number samples, number variants) is small enough, then the standard SKAT p-value calculation is used; (2) if min(number samples, number variants) is too large for standard SKAT, but small enough to explicitly compute the genotype covariance matrix, random projections are used to approximate the eigenvalues of the covariance matrix, and the fastSKAT p-value calculation is used; (3) if min(number samples, number variants) is too big to explicitly compute the genotype covariance matrix, random projections are used to approximate both the eigenvalues and the trace of the covariance matrix, and the fastSKAT p-value calculation is used.)
The fastSMMAT method uses the same random matrix theory as fastSKAT to speed up the computation of the p-value for the adjusted SKAT component of the test. When test = "fastSMMAT"
, the function uses the same logic as for fastSKAT to determine which p-value calculation approach to use for each aggregation unit.
The BinomiRare test, run by using test = "BinomiRare"
, and the CMP test, run by using test = "CMP"
are carrier-only, robust tests. Only variants where the effect allele is minor will be tested. Both tests focuse on carriers of the rare variant allele ("carriers"), and use the estimated probabilities of the binary outcome within the carriers, estimated under the null of not association, and the actual number of observed outcomes, to compute p-values. BinomiRare uses the Poisson-Binomial distribution, and the CMP uses the Conway-Maxwell-Poisson distribution, and is specifically designed for mixed models. (If test = "CMP"
but null.model$family$mixedmodel = FALSE
, the BinomiRare test will be run instead.) These tests provide both a traditional p-value ("pval"
) and a mid-p-value ("midp"
), which is less conservative/more liberal, with the difference being more pronounced for small number of carriers. The BinomiRare test is described in Sofer (2017) and compared to the Score and SPA in various settings in Sofer and Guo (2020).
p-values that are calculated using pchisq
and are smaller than .Machine\$double.xmin
are set to .Machine\$double.xmin
.
A list with the following items:
results |
A data.frame containing the results from the main analysis. Each row is a separate aggregate test: |
If gdsobj
is a SeqVarWindowIterator
:
chr |
The chromosome value |
start |
The start position of the window |
end |
The end position of the window |
Always:
n.site |
The number of variant sites included in the test. |
n.alt |
The number of alternate (effect) alleles included in the test. |
n.sample.alt |
The number of samples with an observed alternate (effect) allele at any variant in the aggregate set. |
If test
is "Burden"
:
Score |
The value of the score function |
Score.SE |
The estimated standard error of the Score |
Score.Stat |
The score Z test statistic |
Score.pval |
The score p-value |
Est |
An approximation of the effect size estimate for each additional unit of burden |
Est.SE |
An approximation of the standard error of the effect size estimate |
PVE |
An approximation of the proportion of phenotype variance explained |
If test
is "SKAT"
or "fastSKAT"
:
Q |
The SKAT test statistic. |
pval |
The SKAT p-value. |
err |
Takes value 1 if there was an error in calculating the p-value; takes the value 2 when multiple random projections were required to get a good approximation from fastSKAT (the reported p-value is likely still reliable); 0 otherwise. |
pval.method |
The p-value calculation method used. When standard SKAT is used, one of "integration" or "saddlepoint"; when fastSKAT random projections are used to approximate eigenvalues of the genotype covariance matrix, one of "ssvd_integration" or "ssvd_saddlepoint"; when fastSKAT random projections are used to approximate both the eigenvalues and the trace of the genotype covariance matrix, one of "rsvd_integration" or "rsvd_saddlepoint". |
If test
is "SMMAT"
or "fastSMMAT"
:
Score_burden |
The value of the score function for the burden test |
Score.SE_burden |
The estimated standard error of the Score for the burden test |
Stat_burden |
The score Z test statistic for the burden test |
pval_burden |
The burden test p-value. |
Q_theta |
The test statistic for the adjusted SKAT test (which is asymptotically independent of the burden test) |
pval_theta |
The p-value of the adjusted SKAT test (which is asymptotically independent of the burden test) |
pval_SMMAT |
The SMMAT p-value after combining pval_burden and pval_theta using Fisher's method. |
err |
Takes value 1 if there was an error calculating the SMMAT p-value; 0 otherwise. If |
pval_theta.method |
The p-value calculation method used for |
If test
is "SKATO"
:
Q_rho |
The SKAT test statistic for the value of rho specified. There will be as many of these variables as there are rho values chosen. |
pval_rho |
The SKAT p-value for the value of rho specified. There will be as many of these variables as there are rho values chosen. |
err_rho |
Takes value 1 if there was an error in calculating the p-value for the value of rho specified when using the "kuonen" or "davies" methods; 0 otherwise. When there is an error, the p-value returned is from the "liu" method. There will be as many of these variables as there are rho values chosen. |
min.pval |
The minimum p-value among the p-values calculated for each choice of rho. |
opt.rho |
The optimal rho value; i.e. the rho value that gave the minimum p-value. |
pval_SKATO |
The SKAT-O p-value after adjustment for searching across multiple rho values. |
If test
is "BinomiRare" or "CMP"
:
n.carrier |
Number of individuals with at least one copy of the effect allele |
n.D.carrier |
Number of cases with at least one copy of the effect allele |
pval |
p-value |
mid.pval |
mid-p-value |
variantInfo |
A list with as many elements as aggregate tests performed. Each element of the list is a data.frame providing information on the variants used in the aggregate test with results presented in the corresponding row of |
variant.id |
The variant ID |
chr |
The chromosome value |
pos |
The base pair position |
allele.index |
The index of the alternate allele. For biallelic variants, this will always be 1. |
n.obs |
The number of samples with non-missing genotypes |
freq |
The estimated effect allele frequency |
MAC |
The minor allele count. For multiallelic variants, "minor" is determined by comparing the count of the allele specified by |
weight |
The weight assigned to the variant in the analysis. |
Matthew P. Conomos, Stephanie M. Gogarten, Thomas Lumley, Tamar Sofer, Ken Rice, Chaoyu Yu, Han Chen
Leal, S.M. & Li, B. (2008). Methods for Detecting Associations with Rare Variants for Common Diseases: Application to Analysis of Sequence Data. American Journal of Human Genetics, 83(3): 311-321.
Browning, S.R. & Madsen, B.E. (2009). A Groupwise Association Test for Rare Mutations Using a Weighted Sum Statistic. PLoS Genetics, 5(2): e1000384.
Wu, M.C, Lee, S., Cai, T., Li, Y., Boehnke, M., & Lin, X. (2011). Rare-Variant Association Testing for Sequencing Data with the Sequence Kernel Association Test. American Journal of Human Genetics, 89(1): 82-93.
Lee, S. et al. (2012). Optimal Unified Approach for Rare-Variant Association Testing with Application to Small-Sample Case-Control Whole-Exome Sequencing Studies. American Journal of Human Genetics, 91(2): 224-237.
Chen, H., Huffman, J. E., Brody, J. A., Wang, C., Lee, S., Li, Z., ... & Blangero, J. (2019). Efficient variant set mixed model association tests for continuous and binary traits in large-scale whole-genome sequencing studies. The American Journal of Human Genetics, 104(2), 260-274.
Lumley, T., Brody, J., Peloso, G., Morrison, A., & Rice, K. (2018). FastSKAT: Sequence kernel association tests for very large sets of markers. Genetic epidemiology, 42(6), 516-527.
library(SeqVarTools)
library(Biobase)
library(GenomicRanges)
# open a sequencing GDS file
gdsfile <- seqExampleFileName("gds")
gds <- seqOpen(gdsfile)
# simulate some phenotype data
set.seed(4)
data(pedigree)
pedigree <- pedigree[match(seqGetData(gds, "sample.id"), pedigree$sample.id),]
pedigree$outcome <- rnorm(nrow(pedigree))
# construct a SeqVarData object
seqData <- SeqVarData(gds, sampleData=AnnotatedDataFrame(pedigree))
# fit the null model
nullmod <- fitNullModel(seqData, outcome="outcome", covars="sex")
# burden test - Range Iterator
gr <- GRanges(seqnames=rep(1,3), ranges=IRanges(start=c(1e6, 2e6, 3e6), width=1e6))
iterator <- SeqVarRangeIterator(seqData, variantRanges=gr)
assoc <- assocTestAggregate(iterator, nullmod, test="Burden",
BPPARAM=BiocParallel::SerialParam())
assoc$results
lapply(assoc$variantInfo, head)
# SKAT test - Window Iterator
seqSetFilterChrom(seqData, include="22")
iterator <- SeqVarWindowIterator(seqData)
assoc <- assocTestAggregate(iterator, nullmod, test="SKAT",
BPPARAM=BiocParallel::SerialParam())
head(assoc$results)
head(assoc$variantInfo)
# SKAT-O test - List Iterator
seqResetFilter(iterator)
gr <- GRangesList(
GRanges(seqnames=rep(22,2), ranges=IRanges(start=c(16e6, 17e6), width=1e6)),
GRanges(seqnames=rep(22,2), ranges=IRanges(start=c(18e6, 20e6), width=1e6)))
iterator <- SeqVarListIterator(seqData, variantRanges=gr)
assoc <- assocTestAggregate(iterator, nullmod, test="SKAT", rho=seq(0, 1, 0.25),
BPPARAM=BiocParallel::SerialParam())
assoc$results
assoc$variantInfo
# user-specified weights - option 1
seqResetFilter(iterator)
variant.id <- seqGetData(gds, "variant.id")
weights <- data.frame(variant.id, weight=runif(length(variant.id)))
variantData(seqData) <- AnnotatedDataFrame(weights)
iterator <- SeqVarListIterator(seqData, variantRanges=gr)
assoc <- assocTestAggregate(iterator, nullmod, test="Burden", weight.user="weight",
BPPARAM=BiocParallel::SerialParam())
assoc$results
assoc$variantInfo
# user-specified weights - option 2
seqResetFilter(iterator)
variantData(seqData)$weight <- NULL
gr <- GRangesList(
GRanges(seqnames=rep(22,2), ranges=IRanges(start=c(16e6, 17e6), width=1e6), weight=runif(2)),
GRanges(seqnames=rep(22,2), ranges=IRanges(start=c(18e6, 20e6), width=1e6), weight=runif(2)))
iterator <- SeqVarListIterator(seqData, variantRanges=gr)
assoc <- assocTestAggregate(iterator, nullmod, test="Burden", weight.user="weight",
BPPARAM=BiocParallel::SerialParam())
assoc$results
assoc$variantInfo
seqClose(seqData)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.