This function calculates the Jensen-Shannon Divergence (JSD) in a
SummarizedExperiment
object.
# S4 method for class 'ANY'
calculateJSD(x, ...)
# S4 method for class 'SummarizedExperiment'
calculateJSD(
x,
assay.type = assay_name,
assay_name = exprs_values,
exprs_values = "counts",
transposed = FALSE,
...
)
runJSD(x, BPPARAM = SerialParam(), chunkSize = nrow(x))
a numeric matrix with samples as rows or a
SummarizedExperiment
object.
optional arguments not used.
Character scalar
. Specifies the name of the
assay used in calculation. (Default: "counts"
)
Deprecated. Use assay.type
instead.
Deprecated. Use assay.type
instead.
Logical scalar
. Is x transposed with samples in rows?
(Default: FALSE
)
A
BiocParallelParam
object specifying whether the calculation should be parallelized.
Integer scalar
. Defines the size of data send
to the individual worker. Only has an effect, if BPPARAM
defines
more than one worker. (Default: nrow(x)
)
a sample-by-sample distance matrix, suitable for NMDS, etc.
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
data(enterotype)
library(scater)
#> Loading required package: scuttle
#> Loading required package: ggplot2
jsd <- calculateJSD(enterotype)
class(jsd)
#> [1] "dist"
head(jsd)
#> [1] 0.1282428 0.1438455 0.1001081 0.2694086 0.2413318 0.2186908
enterotype <- runMDS(enterotype, FUN = calculateJSD, name = "JSD",
exprs_values = "counts")
head(reducedDim(enterotype))
#> [,1] [,2]
#> AM.AD.1 -0.2993825 -0.14425004
#> AM.AD.2 -0.1876141 -0.13490508
#> AM.F10.T1 -0.1761079 -0.02830419
#> AM.F10.T2 -0.1508764 -0.06346108
#> DA.AD.1 -0.2822249 -0.11871798
#> DA.AD.1T -0.3040909 -0.11856633
head(attr(reducedDim(enterotype),"eig"))
#> [1] 12.8663937 2.9088395 1.8794319 0.8494149 0.6543775 0.5578214
attr(reducedDim(enterotype),"GOF")
#> [1] 0.5602642 0.6867101