calculateDMN | R Documentation |
These functions are accessors for functions implemented in the
DirichletMultinomial
package
calculateDMN(x, ...)
## S4 method for signature 'ANY'
calculateDMN(
x,
k = 1,
BPPARAM = SerialParam(),
seed = runif(1, 0, .Machine$integer.max),
...
)
## S4 method for signature 'SummarizedExperiment'
calculateDMN(
x,
assay.type = assay_name,
assay_name = exprs_values,
exprs_values = "counts",
transposed = FALSE,
...
)
runDMN(x, name = "DMN", ...)
getDMN(x, name = "DMN", ...)
## S4 method for signature 'SummarizedExperiment'
getDMN(x, name = "DMN")
bestDMNFit(x, name = "DMN", type = c("laplace", "AIC", "BIC"), ...)
## S4 method for signature 'SummarizedExperiment'
bestDMNFit(x, name = "DMN", type = c("laplace", "AIC", "BIC"))
getBestDMNFit(x, name = "DMN", type = c("laplace", "AIC", "BIC"), ...)
## S4 method for signature 'SummarizedExperiment'
getBestDMNFit(x, name = "DMN", type = c("laplace", "AIC", "BIC"))
calculateDMNgroup(x, ...)
## S4 method for signature 'ANY'
calculateDMNgroup(
x,
variable,
k = 1,
seed = runif(1, 0, .Machine$integer.max),
...
)
## S4 method for signature 'SummarizedExperiment'
calculateDMNgroup(
x,
variable,
assay.type = assay_name,
assay_name = exprs_values,
exprs_values = "counts",
transposed = FALSE,
...
)
performDMNgroupCV(x, ...)
## S4 method for signature 'ANY'
performDMNgroupCV(
x,
variable,
k = 1,
seed = runif(1, 0, .Machine$integer.max),
...
)
## S4 method for signature 'SummarizedExperiment'
performDMNgroupCV(
x,
variable,
assay.type = assay_name,
assay_name = exprs_values,
exprs_values = "counts",
transposed = FALSE,
...
)
x |
a numeric matrix with samples as rows or a
|
... |
optional arguments not used. |
k |
|
BPPARAM |
A
|
seed |
|
assay.type |
|
assay_name |
Deprecated. Use |
exprs_values |
Deprecated. Use |
transposed |
|
name |
|
type |
|
variable |
|
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
# Get and display which objects fits best
bestFit <- metadata(tse)$DMM$best
bestFit
# Get the model that generated the best fit
bestModel <- metadata(tse)$DMM$dmm[[bestFit]]
bestModel
# Get the sample-cluster assignment probability matrix
head(metadata(tse)$DMM$prob)
# Get the weight of each component for the best model
bestModel@mixture$Weight
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.