R/inference_methods_generic.R
dcScore.Rd
Implementations and wrappers for existing implementations for methods inferring differential associations/co-expression. This method requires a matrix of expression and a binary condition to compute the differential association scores for all pairs of features (genes). Applications are not limited to analysis of gene expression data and may be used for differential associations in general.
dcScore(emat, condition, dc.method, ...)
# S4 method for matrix
dcScore(emat, condition, dc.method = "zscore", ...)
# S4 method for Matrix
dcScore(emat, condition, dc.method = "zscore", ...)
# S4 method for data.frame
dcScore(emat, condition, dc.method = "zscore", ...)
# S4 method for ExpressionSet
dcScore(emat, condition, dc.method = "zscore", ...)
# S4 method for SummarizedExperiment
dcScore(emat, condition, dc.method = "zscore", ...)
# S4 method for DGEList
dcScore(emat, condition, dc.method = "zscore", ...)
a matrix, Matrix, data.frame, ExpressionSet, SummarizedExperiment or DGEList
a numeric, (with 1's and 2's representing a binary condition), a factor with 2 levels or a character representing 2 conditions
a character, representing the method to use. Use
dcMethods()
to get a list of methods
possible arguments are cor.method
, diffcoex.beta
,
ebcoexpress.useBWMC
, ebcoexpress.plot
, ldgm.lambda
,
ldgm.ntarget
and ldgm.iter
. See details
a matrix, of scores/statistics representing differential associations; p-values will be returned if FTGI is used and posterior probabilities if EBcoexpress is used.
When using data from sequencing experiments, make sure appropriate filtering for low counts and data transformation has been performed. Not doing so will affect estimation of correlation coefficients which most methods rely on.
Additional method specific parameters can be supplied to the function.
cor.method
can be set to either 'pearson' (default) or 'spearman' to
determine the method to use for estimating correlations. These are the two
measures currently supported in the package. We recommend using the
'spearman' correlation when dealing with sequencing data.
The beta parameter in the DiffCoEx method can be specified using
diffcoex.beta
(defaults to 6). This enable soft thresholding of
correlations similar to WGCNA.
EBcoexpress specific parameters include ebcoexpress.useBWMC
(defaults to TRUE
) representing whether to use the bi-weight
mid-correlation coefficient or not, and ebcoexpress.plot
which plots
the diagnostic plots if set to TRUE
(defaults to FALSE
).
LDGM specific parameters include ldgm.lambda
, ldgm.ntarget
and ldgm.iter
. ldgm.lambda
specifies the L1 regularisation
parameter to use when fitting the model. This can be tuned and specified by
the user. Alternatively, this can be tuned such that the resulting network
has a specified number of edges. In this case, ldgm.ntarget
should
be specified instead. ldgm.iter
is the maximum number of iterations
to perform when tuning ldgm.lambda
using ldgm.ntarget
(defaults to 50).
EBcoexpress and GGM-based are implemented by providing interfaces to,
or using functions from the EBcoexpress
, GeneNet
, and
COSINE
packages respectively. If using any of these methods, please
cite the appropriate packages and the appropriate methodology
articles.
x <- matrix(rnorm(60), 2, 30)
cond <- rep(1:2, 15)
dcScore(x, cond) #defaults to zscore
#> 1 2
#> 1 NA -0.5583278
#> 2 -0.5583278 NA
#> attr(,"cor.method")
#> [1] "pearson"
#> attr(,"call")
#> z.score(emat = emat, condition = condition)
#> attr(,"dc.method")
#> [1] "zscore"
dcScore(x, cond, dc.method = 'diffcoex')
#> 1 2
#> 1 NA 0.000431787
#> 2 0.000431787 NA
#> attr(,"cor.method")
#> [1] "pearson"
#> attr(,"diffcoex.beta")
#> [1] 6
#> attr(,"call")
#> diffcoex.score(emat = emat, condition = condition)
#> attr(,"dc.method")
#> [1] "diffcoex"