getCrossAssociation | R Documentation |
Calculate correlations between features of two experiments.
getCrossAssociation(x, ...)
## S4 method for signature 'MultiAssayExperiment'
getCrossAssociation(
x,
experiment1 = 1,
experiment2 = 2,
assay.type1 = assay_name1,
assay_name1 = "counts",
assay.type2 = assay_name2,
assay_name2 = "counts",
altexp1 = NULL,
altexp2 = NULL,
col.var1 = colData_variable1,
colData_variable1 = NULL,
col.var2 = colData_variable2,
colData_variable2 = NULL,
by = MARGIN,
MARGIN = 1,
method = c("kendall", "spearman", "categorical", "pearson"),
mode = "table",
p.adj.method = p_adj_method,
p_adj_method = c("fdr", "BH", "bonferroni", "BY", "hochberg", "holm", "hommel", "none"),
p.adj.threshold = p_adj_threshold,
p_adj_threshold = NULL,
cor.threshold = cor_threshold,
cor_threshold = NULL,
sort = FALSE,
filter.self.cor = filter_self_correlations,
filter_self_correlations = FALSE,
verbose = TRUE,
test.signif = test_significance,
test_significance = FALSE,
show.warnings = show_warnings,
show_warnings = TRUE,
paired = FALSE,
...
)
## S4 method for signature 'SummarizedExperiment'
getCrossAssociation(x, experiment2 = x, ...)
x |
A
|
... |
Additional arguments:
|
experiment1 |
|
experiment2 |
|
assay.type1 |
|
assay_name1 |
Deprecated. Use |
assay.type2 |
|
assay_name2 |
Deprecated. Use |
altexp1 |
|
altexp2 |
|
col.var1 |
|
colData_variable1 |
Deprecated. Use |
col.var2 |
|
colData_variable2 |
Deprecated. Use |
by |
A |
MARGIN |
Deprecated. Use |
method |
|
mode |
|
p.adj.method |
|
p_adj_method |
Deprecated. Use |
p.adj.threshold |
|
p_adj_threshold |
Deprecated. Use |
cor.threshold |
|
cor_threshold |
Deprecated. Use |
sort |
|
filter.self.cor |
|
filter_self_correlations |
Deprecated. Use |
verbose |
|
test.signif |
|
test_significance |
Deprecated. Use |
show.warnings |
|
show_warnings |
Deprecated. use |
paired |
|
The function getCrossAssociation
calculates associations between
features of two experiments. By default, it not only computes associations
but also tests their significance. If desired, setting
test.signif
to FALSE disables significance calculation.
We recommend the non-parametric Kendall's tau as the default method for association analysis. Kendall's tau has desirable statistical properties and robustness at lower sample sizes. Spearman rank correlation can provide faster solutions when running times are critical.
This function returns associations in table or matrix format. In table format, returned value is a data frame that includes features and associations (and p-values) in columns. In matrix format, returned value is a one matrix when only associations are calculated. If also significances are tested, then returned value is a list of matrices.
data(HintikkaXOData)
mae <- HintikkaXOData
# Subset so that less observations / quicker to run, just for example
mae[[1]] <- mae[[1]][1:20, 1:10]
mae[[2]] <- mae[[2]][1:20, 1:10]
# Several rows in the counts assay have a standard deviation of zero
# Remove them, since they do not add useful information about
# cross-association
mae[[1]] <- mae[[1]][rowSds(assay(mae[[1]])) > 0, ]
# Transform data
mae[[1]] <- transformAssay(mae[[1]], method = "rclr")
# Calculate cross-correlations
result <- getCrossAssociation(mae, method = "pearson", assay.type2 = "nmr")
# Show first 5 entries
head(result, 5)
# Use altExp option to specify alternative experiment from the experiment
altExp(mae[[1]], "Phylum") <- agglomerateByRank(mae[[1]], rank = "Phylum")
# Transform data
altExp(mae[[1]], "Phylum") <- transformAssay(
altExp(mae[[1]], "Phylum"), method = "relabundance")
# When mode = "matrix", the return value is a matrix
result <- getCrossAssociation(
mae, experiment2 = 2, assay.type1 = "relabundance", assay.type2 = "nmr",
altexp1 = "Phylum", method = "pearson", mode = "matrix")
# Show first 5 entries
head(result, 5)
# If test.signif = TRUE, then getCrossAssociation additionally returns
# significances
# filter.self.cor = TRUE filters self correlations
# p.adj.threshold can be used to filter those features that do not
# have any correlations whose p-value is lower than the threshold
result <- getCrossAssociation(
mae[[1]], experiment2 = mae[[1]], method = "pearson",
filter.self.cor = TRUE, p.adj.threshold = 0.05, test.signif = TRUE)
# Show first 5 entries
head(result, 5)
# Returned value is a list of matrices
names(result)
# Calculate Bray-Curtis dissimilarity between samples. If dataset includes
# paired samples, you can use paired = TRUE.
result <- getCrossAssociation(
mae[[1]], mae[[1]], by = 2, paired = FALSE,
association.fun = getDissimilarity, method = "bray")
# If experiments are equal and measure is symmetric
# (e.g., taxa1 vs taxa2 == taxa2 vs taxa1),
# it is possible to speed-up calculations by calculating association only
# for unique variable-pairs. Use "symmetric" to choose whether to measure
# association for only other half of of variable-pairs.
result <- getCrossAssociation(
mae, experiment1 = "microbiota", experiment2 = "microbiota",
assay.type1 = "counts", assay.type2 = "counts", symmetric = TRUE)
# For big data sets, the calculations might take a long time.
# To speed them up, you can take a random sample from the data.
# When dealing with complex biological problems, random samples can be
# enough to describe the data. Here, our random sample is 30 % of whole data.
sample_size <- 0.3
tse <- mae[[1]]
tse_sub <- tse[ sample( seq_len( nrow(tse) ), sample_size * nrow(tse) ), ]
result <- getCrossAssociation(tse_sub)
# It is also possible to choose variables from colData and calculate
# association between assay and sample metadata or between variables of
# sample metadata
mae[[1]] <- addAlpha(mae[[1]])
# col.var works similarly to assay.type. Instead of fetching
# an assay named assay.type from assay slot, it fetches a column named
# col.var from colData.
result <- getCrossAssociation(
mae[[1]], assay.type1 = "counts",
col.var2 = c("shannon_diversity", "coverage_diversity"),
test.signif = TRUE)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.