R/AllGenerics.R
, R/calculateDMM.R
calculateDMN.Rd
These functions are accessors for functions implemented in the
DirichletMultinomial
package
calculateDMN(x, ...)
getDMN(x, name = "DMN", ...)
bestDMNFit(x, name = "DMN", type = c("laplace", "AIC", "BIC"), ...)
calculateDMNgroup(x, ...)
performDMNgroupCV(x, ...)
getBestDMNFit(x, name = "DMN", type = c("laplace", "AIC", "BIC"), ...)
# S4 method for class 'ANY'
calculateDMN(
x,
k = 1,
BPPARAM = SerialParam(),
seed = runif(1, 0, .Machine$integer.max),
...
)
# S4 method for class 'SummarizedExperiment'
calculateDMN(
x,
assay.type = assay_name,
assay_name = exprs_values,
exprs_values = "counts",
transposed = FALSE,
...
)
runDMN(x, name = "DMN", ...)
# S4 method for class 'SummarizedExperiment'
getDMN(x, name = "DMN")
# S4 method for class 'SummarizedExperiment'
bestDMNFit(x, name = "DMN", type = c("laplace", "AIC", "BIC"))
# S4 method for class 'SummarizedExperiment'
getBestDMNFit(x, name = "DMN", type = c("laplace", "AIC", "BIC"))
# S4 method for class 'ANY'
calculateDMNgroup(
x,
variable,
k = 1,
seed = runif(1, 0, .Machine$integer.max),
...
)
# S4 method for class 'SummarizedExperiment'
calculateDMNgroup(
x,
variable,
assay.type = assay_name,
assay_name = exprs_values,
exprs_values = "counts",
transposed = FALSE,
...
)
# S4 method for class 'ANY'
performDMNgroupCV(
x,
variable,
k = 1,
seed = runif(1, 0, .Machine$integer.max),
...
)
# S4 method for class 'SummarizedExperiment'
performDMNgroupCV(
x,
variable,
assay.type = assay_name,
assay_name = exprs_values,
exprs_values = "counts",
transposed = FALSE,
...
)
a numeric matrix with samples as rows or a
SummarizedExperiment
object.
optional arguments not used.
Character scalar
. The name to store the result in
metadata
Character scalar
. The type of measure used for the
goodness of fit. One of ‘laplace’, ‘AIC’ or ‘BIC’.
Numeric scalar
. The number of Dirichlet components to fit.
See dmn
. (Default: 1
)
A
BiocParallelParam
object specifying whether the calculation should be parallelized.
Numeric scalar
. Random number seed. See
dmn
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
)
Character scalar
. A variable from colData
to
use as a grouping variable. Must be a character of factor.
calculateDMN
and getDMN
return a list of DMN
objects,
one element for each value of k provided.
bestDMNFit
returns the index for the best fit and getBestDMNFit
returns a single DMN
object.
calculateDMNgroup
returns a
DMNGroup
object
performDMNgroupCV
returns a data.frame
DMN-class
,
DMNGroup-class
,
dmn
,
dmngroup
,
cvdmngroup
,
accessors for DMN objects
fl <- system.file(package="DirichletMultinomial", "extdata", "Twins.csv")
counts <- as.matrix(read.csv(fl, row.names=1))
fl <- system.file(package="DirichletMultinomial", "extdata", "TwinStudy.t")
pheno0 <- scan(fl)
lvls <- c("Lean", "Obese", "Overwt")
pheno <- factor(lvls[pheno0 + 1], levels=lvls)
colData <- DataFrame(pheno = pheno)
tse <- TreeSummarizedExperiment(assays = list(counts = counts),
colData = colData)
library(bluster)
# Compute DMM algorithm and store result in metadata
tse <- addCluster(tse, name = "DMM", DmmParam(k = 1:3, type = "laplace"),
by = "samples", full = TRUE)
# Get the list of DMN objects
metadata(tse)$DMM$dmm
#> [[1]]
#> class: DMN
#> k: 1
#> samples x taxa: 278 x 130
#> Laplace: 39227.38 BIC: 39527.91 AIC: 39292.11
#>
#> [[2]]
#> class: DMN
#> k: 2
#> samples x taxa: 278 x 130
#> Laplace: 38872.71 BIC: 39588.93 AIC: 39115.53
#>
#> [[3]]
#> class: DMN
#> k: 3
#> samples x taxa: 278 x 130
#> Laplace: 38854.24 BIC: 40001.15 AIC: 39290.14
#>
# Get and display which objects fits best
bestFit <- metadata(tse)$DMM$best
bestFit
#> [1] 3
# Get the model that generated the best fit
bestModel <- metadata(tse)$DMM$dmm[[bestFit]]
bestModel
#> class: DMN
#> k: 3
#> samples x taxa: 278 x 130
#> Laplace: 38854.24 BIC: 40001.15 AIC: 39290.14
# Get the sample-cluster assignment probability matrix
head(metadata(tse)$DMM$prob)
#> 1 2 3
#> TS1.2 9.999997e-01 2.999201e-07 4.862193e-09
#> TS10.2 9.516058e-01 8.378014e-09 4.839415e-02
#> TS100.2 4.436106e-08 8.567014e-02 9.143298e-01
#> TS100 9.919489e-01 8.050638e-03 4.122346e-07
#> TS101.2 1.352493e-12 4.164390e-07 9.999996e-01
#> TS103.2 9.999977e-01 7.471442e-10 2.313071e-06
# Get the weight of each component for the best model
bestModel@mixture$Weight
#> [1] 0.5655825 0.2227240 0.2116935