addDissimilarity | R Documentation |
These functions are designed to calculate dissimilarities on data stored
within a
TreeSummarizedExperiment
object. For overlap, Unifrac, and Jensen-Shannon Divergence (JSD)
dissimilarities, the functions use mia internal functions, while for other
types of dissimilarities, they rely on vegdist
by default.
addDissimilarity(x, method, ...)
getDissimilarity(x, method, ...)
## S4 method for signature 'SummarizedExperiment'
addDissimilarity(x, method = "bray", name = method, ...)
## S4 method for signature 'SummarizedExperiment'
getDissimilarity(
x,
method = "bray",
assay.type = "counts",
niter = NULL,
transposed = FALSE,
...
)
## S4 method for signature 'TreeSummarizedExperiment'
getDissimilarity(
x,
method = "bray",
assay.type = "counts",
niter = NULL,
transposed = FALSE,
...
)
## S4 method for signature 'ANY'
getDissimilarity(x, method = "bray", niter = NULL, ...)
x |
|
method |
|
... |
other arguments passed into
|
name |
|
assay.type |
|
niter |
|
transposed |
|
Overlap reflects similarity between sample-pairs. When overlap is calculated using relative abundances, the higher the value the higher the similarity is. When using relative abundances, overlap value 1 means that all the abundances of features are equal between two samples, and 0 means that samples have completely different relative abundances.
Unifrac is calculated with rbiom:unifrac()
.
If rarefaction is enabled, vegan:avgdist()
is
utilized.
Rarefaction can be used to control uneven sequencing depths. Although, it is highly debated method. Some think that it is the only option that successfully controls the variation caused by uneven sampling depths. The biggest argument against rarefaction is the fact that it omits data.
Rarefaction works by sampling the counts randomly. This random sampling
is done niter
times. In each sampling iteration, sample
number
of random samples are drawn, and dissimilarity is calculated for this
subset. After the iterative process, there are niter
number of
result that are then averaged to get the final result.
Refer to Schloss (2024) for more details on rarefaction.
getDissimilarity
returns a sample-by-sample dissimilarity matrix.
addDissimilarity
returns x
that includes dissimilarity matrix
in its metadata.
For unifrac dissimilarity: http://bmf.colorado.edu/unifrac/
See also additional descriptions of Unifrac in the following articles:
Lozupone, Hamady and Knight, “Unifrac - An Online Tool for Comparing Microbial Community Diversity in a Phylogenetic Context.”, BMC Bioinformatics 2006, 7:371
Lozupone, Hamady, Kelley and Knight, “Quantitative and qualitative (beta) diversity measures lead to different insights into factors that structure microbial communities.” Appl Environ Microbiol. 2007
Lozupone C, Knight R. “Unifrac: a new phylogenetic method for comparing microbial communities.” Appl Environ Microbiol. 2005 71 (12):8228-35.
For JSD dissimilarity: Jensen-Shannon Divergence and Hilbert space embedding. Bent Fuglede and Flemming Topsoe University of Copenhagen, Department of Mathematics http://www.math.ku.dk/~topsoe/ISIT2004JSD.pdf
For rarefaction: Schloss PD (2024) Rarefaction is currently the best approach to control for uneven sequencing effort in amplicon sequence analyses. mSphere 28;9(2):e0035423. doi: 10.1128/msphere.00354-23
http://en.wikipedia.org/wiki/Jensen-Shannon_divergence
library(mia)
library(scater)
# load dataset
data(GlobalPatterns)
tse <- GlobalPatterns
### Overlap dissimilarity
tse <- addDissimilarity(tse, method = "overlap", detection = 0.25)
metadata(tse)[["overlap"]][1:6, 1:6]
### JSD dissimilarity
tse <- addDissimilarity(tse, method = "jsd")
metadata(tse)[["jsd"]][1:6, 1:6]
# Multi Dimensional Scaling applied to JSD dissimilarity matrix
tse <- addMDS(tse, method = "overlap", assay.type = "counts")
reducedDim(tse, "MDS") |> head()
### Unifrac dissimilarity
res <- getDissimilarity(tse, method = "unifrac", weighted = FALSE)
dim(as.matrix(res))
tse <- addDissimilarity(tse, method = "unifrac", weighted = TRUE)
metadata(tse)[["unifrac"]][1:6, 1:6]
### Bray dissimilarity
# Bray is usually applied to relative abundances so we have to apply
# transformation first
tse <- transformAssay(tse, method = "relabundance")
res <- getDissimilarity(tse, method = "bray", assay.type = "relabundance")
as.matrix(res)[1:6, 1:6]
# If applying rarefaction, the input must be count matrix and transformation
# method specified in function call (Note: increase niter)
rclr <- function(x){
vegan::decostand(x, method="rclr")
}
res <- getDissimilarity(
tse, method = "euclidean", transf = rclr, niter = 2L)
as.matrix(res)[1:6, 1:6]
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.