compare_conditions | R Documentation |
Take modules built from multiples conditions and search for preservation, non-preservation or one of them, against one or mutliple conditions of reference. Use 7 topological features to perform the differents test, and use permutation to validate results.
compare_conditions(
data_expr_list,
adja_list,
cor_list = NULL,
modules_list,
ref = names(data_expr_list)[1],
test = NULL,
cor_func = c("pearson", "spearman", "bicor", "other"),
your_func = NULL,
n_perm = 10000,
test_alternative_hyp = c("greater", "less", "two.sided"),
pvalue_th = 0.01,
n_threads = NULL,
...
)
data_expr_list |
list of matrix or data.frame or SummarizedExperiment, list of expression data by condition, with genes as column and samples as row. |
adja_list |
list of adjacency matrices, list of square tables by condition, representing connectivity between each genes as returned by build_net. |
cor_list |
list of matrices and/or data.frames, list of square tables
by condition, representing correlation between each gene. Must be the same
used to create networks in |
modules_list |
list of modules or nested list of modules, list of modules in one condition (will be considered as the one from reference) or a condition named list with list of modules built in each one. |
ref |
string or vector of strings, condition(s) name to be used as reference for permutation tests, or "cross comparison" if you want to compare each condition with the other as reference. Default will be the name of the first element in data_expr_list. |
test |
string or vector of strings, condition(s) name to be tested for permutation tests. If NULL, all conditions except these in ref will be taken. If ref is set to "cross comparison", any test specified will be ignored. |
cor_func |
string, name of the correlation function to be used. Must be one of "pearson", "spearman", "bicor", "other". If "other", your_func must be provided |
your_func |
function returning correlation values. Final values must be in [-1;1] |
n_perm |
integer, number of permutation, meaning number of random gene name re-assignment inside network to compute all tests and statistics for module comparison between condition. |
test_alternative_hyp |
string, either "greater", "less" or "two.sided".
Alternative hypothesis (H1) used for the permutation test. Determine if the
metrics computed on permuted values are expected to be greater, less or both
than the observed ones. More details: |
pvalue_th |
decimal, threshold of pvalue below which test_alternative_hyp is considered significant. If "two.sided", then pvalue_th is splitted in two for each side (preserved/not preserved). |
n_threads |
integer, number of threads that can be used to paralellise the computing |
... |
any other parameter compatible with
|
keep_cor_mat
in build_net
to TRUE.modulePreservation
.A nested list where first element is each ref provided, second level each condition to test, and then elements containing information on the comparison. See NetRep::modulePreservation() for more detail.
expr_by_cond <- list(cond1 = kuehne_expr[1:24, 1:350],
cond2 = kuehne_expr[25:48, 1:350])
net_by_cond <- lapply(expr_by_cond, build_net, cor_func = "spearman",
n_threads = 1, keep_matrices = "both")
mod_by_cond <- mapply(detect_modules, expr_by_cond,
lapply(net_by_cond, `[[`, "network"),
MoreArgs = list(detailled_result = TRUE),
SIMPLIFY = FALSE)
comparison <- compare_conditions(expr_by_cond,
lapply(net_by_cond, `[[`, "adja_mat"),
lapply(net_by_cond, `[[`, "cor_mat"),
lapply(mod_by_cond, `[[`, "modules"),
n_perm = 100)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.