GLMMMiRKAT: The microbiome regression-based kernel association test based...

Description Usage Arguments Value Author(s) References Examples

Description

This function tests the association betweem microbial community (e.g., bacterial kingdom) composition and a host trait of interest (e.g., health/disease status, medical intervention, behavioral/environmental factor) based on correlated (e.g., family-based or longitudinal) microbiome studies. For the host trait of interest, Gaussian (e.g., body mass index), Binomial (e.g., disease status, treatment/placebo) or Poisson (e.g., number of tumors/treatments) traits can be handled.

Usage

1
GLMMMiRKAT(y, cov = NULL, id, time.pt = NULL, Ks, model, slope = FALSE, n.perm = 5000)

Arguments

y

A numeric vector of Gaussian (e.g., body mass index), Binomial (e.g., disease status, treatment/placebo) or Poisson (e.g., number of tumors/treatments) traits.

cov

A data.frame for covariate (e.g., age, gender) adjustment(s). Default is cov = NULL for no covariate adjustment.

id

A vector of cluster (e.g., family, subject including repeated measurements) IDs.

time.pt

A vector of time points for the longitudinal studies. 'time.pt' is not required (i.e., 'time.pt = NULL') for the random intercept model. Default is time.pt = NULL.

Ks

A list of the kernel matrices created using GLMMMiRKAT::Kernels (?Kernels).

model

"gaussian" for Gaussian outcomes, "binomial" for Binomial outcomes, "poisson" for Poisson outcomes.

slope

An indicator to include random slopes in the model (slope = TRUE) or not (slope = FALSE). 'slope = FALSE' is for the random intercept model. 'slope = TRUE' is for the random slope model. For the random slope model (slope = TRUE), 'time.pt' is required.

n.perm

A number of permutations. Default is 5000.

Value

ItembyItem: A vector of the estimated p-values for the item-by-item GLMM-MiRKAT analyses

aGLMMMiRKAT: The estimated p-value for the adaptive GLMM-MiRKAT (aGLMM-MiRKAT) analysis

Author(s)

Hyunwook Koh

References

Koh H, Li Y, Zhan X, Chen J, Zhao N. (2019) A distance-based kernel association test based on the generalized linear mixed model for correlated microbiome studies. Front. Genet. 458(10), 1-14.

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
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
# Import requisite R packages
require(CompQuadForm)
require(dirmult)
require(ecodist)
require(GUniFrac)
require(lme4)
require(MASS)
require(Matrix)
require(permute)
require(phyloseq)

## Example 1. Gaussian (e.g., body mass index) traits
		
# Import example microbiome data with Gaussian traits
data(nor.phy)

# Rarefy the microbiome data using phyloseq::rarefy_even_depth 
# to control differing total reads per sample (recommended).
set.seed(100)
nor.phy <- rarefy_even_depth(nor.phy, rngseed=TRUE)
		
# Extract microbiome and meta information
otu.tab <- otu_table(nor.phy)
tree <- phy_tree(nor.phy)
meta <- sample_data(nor.phy)
y <- meta$y
id <- meta$id
x1 <- meta$x1
x2 <- meta$x2
covs <- as.data.frame(cbind(x1,x2))
covs[,2] <- as.factor(covs[,2]) 
# 'as.factor()' is needed for categorical covariates.

# Create kernel matrices
Ks <- Kernels(otu.tab, tree)
		
# Run GLMM-MiRKAT
GLMMMiRKAT(y, cov=covs, id=id, Ks=Ks, model="gaussian")
	
## Example 2. Binomial (e.g., disease status, treatment/placebo) traits
		
# Import microbiome data with binomial traits
data(bin.phy)

# Rarefy the microbiome data using phyloseq::rarefy_even_depth 
# to control differing total reads per sample (recommended).
set.seed(100)
bin.phy <- rarefy_even_depth(bin.phy, rngseed=TRUE)

# Extract microbiome and meta information
otu.tab <- otu_table(bin.phy)
tree <- phy_tree(bin.phy)
meta <- sample_data(bin.phy)
y <- meta$y
id <- meta$id
x1 <- meta$x1
x2 <- meta$x2
covs <- as.data.frame(cbind(x1,x2))
covs[,2] <- as.factor(covs[,2]) 
# 'as.factor()' is needed for categorical covariates.
		
# Create kernel matrices
Ks <- Kernels(otu.tab, tree)
		
# Run GLMM-MiRKAT
GLMMMiRKAT(y, cov=covs, id=id, Ks=Ks, model="binomial")

## Example 3. Poisson (e.g., number of tumors/treatments) traits
		
# Import microbiome data with Poisson traits
data(pos.phy)

# Rarefy the microbiome data using phyloseq::rarefy_even_depth 
# to control differing total reads per sample (recommended).
set.seed(100)
pos.phy <- rarefy_even_depth(pos.phy, rngseed=TRUE)

# Extract microbiome and meta information
otu.tab <- otu_table(pos.phy)
tree <- phy_tree(pos.phy)
meta <- sample_data(pos.phy)
y <- meta$y
id <- meta$id
x1 <- meta$x1
x2 <- meta$x2
covs <- as.data.frame(cbind(x1,x2))
covs[,2] <- as.factor(covs[,2]) 
# 'as.factor()' is needed for categorical covariates.

# Create kernel matrices
Ks <- Kernels(otu.tab, tree)
	
# Run GLMM-MiRKAT
GLMMMiRKAT(y, cov=covs, id=id, Ks=Ks, model="poisson")

hk1785/GLMM-MiRKAT documentation built on Feb. 29, 2020, 6:06 p.m.