getExperimentCrossAssociation | R Documentation |
Calculate correlations between features of two experiments.
getExperimentCrossAssociation(x, ...)
## S4 method for signature 'MultiAssayExperiment'
getExperimentCrossAssociation(
x,
experiment1 = 1,
experiment2 = 2,
assay.type1 = assay_name1,
assay_name1 = "counts",
assay.type2 = assay_name2,
assay_name2 = "counts",
altexp1 = NULL,
altexp2 = NULL,
colData_variable1 = NULL,
colData_variable2 = NULL,
MARGIN = 1,
method = c("kendall", "spearman", "categorical", "pearson"),
mode = "table",
p_adj_method = c("fdr", "BH", "bonferroni", "BY", "hochberg", "holm", "hommel", "none"),
p_adj_threshold = NULL,
cor_threshold = NULL,
sort = FALSE,
filter_self_correlations = FALSE,
verbose = TRUE,
test_significance = FALSE,
show_warnings = TRUE,
paired = FALSE,
...
)
## S4 method for signature 'SummarizedExperiment'
getExperimentCrossAssociation(x, experiment2 = x, ...)
testExperimentCrossAssociation(x, ...)
## S4 method for signature 'ANY'
testExperimentCrossAssociation(x, ...)
testExperimentCrossCorrelation(x, ...)
## S4 method for signature 'ANY'
testExperimentCrossCorrelation(x, ...)
getExperimentCrossCorrelation(x, ...)
## S4 method for signature 'ANY'
getExperimentCrossCorrelation(x, ...)
x |
A
|
... |
Additional arguments:
|
experiment1 |
A single character or numeric value for selecting the experiment 1
from |
experiment2 |
A single character or numeric value for selecting the experiment 2
from |
assay.type1 |
A single character value for selecting the
|
assay_name1 |
a single |
assay.type2 |
A single character value for selecting the
|
assay_name2 |
a single |
altexp1 |
A single numeric or character value specifying alternative experiment
from the altExp of experiment 1. If NULL, then the experiment is itself
and altExp option is disabled.
(By default: |
altexp2 |
A single numeric or character value specifying alternative experiment
from the altExp of experiment 2. If NULL, then the experiment is itself
and altExp option is disabled.
(By default: |
colData_variable1 |
A character value specifying column(s) from colData
of experiment 1. If colData_variable1 is used, assay.type1 is disabled.
(By default: |
colData_variable2 |
A character value specifying column(s) from colData
of experiment 2. If colData_variable2 is used, assay.type2 is disabled.
(By default: |
MARGIN |
A single numeric value for selecting if association are calculated
row-wise / for features (1) or column-wise / for samples (2). Must be |
method |
A single character value for selecting association method
('kendall', pearson', or 'spearman' for continuous/numeric; 'categorical' for discrete)
(By default: |
mode |
A single character value for selecting output format
Available formats are 'table' and 'matrix'. (By default: |
p_adj_method |
A single character value for selecting adjustment method of
p-values. Passed to |
p_adj_threshold |
A single numeric value (from 0 to 1) for selecting
adjusted p-value threshold for filtering.
(By default: |
cor_threshold |
A single numeric absolute value (from 0 to 1) for selecting
correlation threshold for filtering.
(By default: |
sort |
A single boolean value for selecting whether to sort features or not
in result matrices. Used method is hierarchical clustering.
(By default: |
filter_self_correlations |
A single boolean value for selecting whether to
filter out correlations between identical items. Applies only when correlation
between experiment itself is tested, i.e., when assays are identical.
(By default: |
verbose |
A single boolean value for selecting whether to get messages about progress of calculation. |
test_significance |
A single boolean value for selecting whether to test statistical significance of associations. |
show_warnings |
A single boolean value for selecting whether to show warnings that might occur when correlations and p-values are calculated. |
paired |
A single boolean value for specifying if samples are paired or not.
|
These functions calculates associations between features of two experiments.
getExperimentCrossAssociation
calculates only associations by default.
testExperimentCrossAssociation
calculates also significance of
associations.
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.
These functions return 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.
Leo Lahti and Tuomas Borman. Contact: microbiome.github.io
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 <- getExperimentCrossAssociation(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 <- getExperimentCrossAssociation(mae, experiment2 = 2,
assay.type1 = "relabundance", assay.type2 = "nmr",
altexp1 = "Phylum",
method = "pearson", mode = "matrix")
# Show first 5 entries
head(result, 5)
# testExperimentCorrelation additionally returns significances
# filter_self_correlations = 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 <- testExperimentCrossAssociation(mae[[1]], experiment2 = mae[[1]], method = "pearson",
filter_self_correlations = TRUE,
p_adj_threshold = 0.05)
# Show first 5 entries
head(result, 5)
# getExperimentCrossAssociation also returns significances when
# test_significance = TRUE
# Warnings can be suppressed by using show_warnings = FALSE
result <- getExperimentCrossAssociation(mae[[1]], experiment2 = mae[[2]], method = "pearson",
assay.type2 = "nmr",
mode = "matrix", test_significance = TRUE,
show_warnings = FALSE)
# 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 <- getExperimentCrossAssociation(mae[[1]], mae[[1]], MARGIN = 2, paired = FALSE,
association_FUN = vegan::vegdist, 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 <- getExperimentCrossAssociation(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 <- testExperimentCrossAssociation(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]] <- estimateDiversity(mae[[1]])
# colData_variable works similarly to assay.type. Instead of fetching an assay
# named assay.type from assay slot, it fetches a column named colData_variable
# from colData.
result <- getExperimentCrossAssociation(mae[[1]], assay.type1 = "counts",
colData_variable2 = c("shannon", "coverage"))
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.