Plot evolutionary signatures across transcriptomes, for multiple expression sets
Source:R/PlotSignatureMultiple.R
PlotSignatureMultiple.Rd
PlotSignatureMultiple is used to compare multiple transcriptomic index profiles (e.g. TAI) over a common process of interest (e.g. a developmental process). The main use case is visualising how removing different subsets of genes from the expression set can perturb the transcriptome index signal.
Usage
PlotSignatureMultiple(
ExpressionSets,
set.labels,
measure = "TAI",
TestStatistic = "FlatLineTest",
modules = NULL,
permutations = 1000,
p.value = TRUE,
shaded.area = FALSE,
xlab = "Ontogeny",
ylab = "Transcriptome Index",
main = "",
legend.title = "Expression Sets",
lwd = 4,
alpha = 0.1,
y.ticks = 3
)
Arguments
- ExpressionSets
a list of PhyloExpressionSet, DivergenceExpressionSet or PolymorphismsExpressionSet objects.
- set.labels
a character vector of labels, one for each given expression set
- measure
type of transcriptome index that shall be computed. E.g.
measure = "TAI"
(Transcriptome Age Index)measure = "TDI"
(Transcriptome Divergence Index)measure = "TPI"
(Transcriptome Polymorphism Index)
- TestStatistic
a string defining the type of test statistics to be used to quantify the statistical significance the present phylotranscriptomics pattern. Possible values can be:
TestStatistic
="FlatLineTest"
: Statistical test for the deviation from a flat lineTestStatistic
="ReductiveHourglassTest"
: Statistical test for the existence of a hourglass shape (high-low-high pattern)TestStatistic
="ReverseHourglassTest"
: Statistical test for the existence of a reverse hourglass pattern (low-high-low pattern)TestStatistic
="EarlyConservationTest"
: Statistical test for the existence of a early conservation pattern (low-high-high pattern)TestStatistic
="LateConservationTest"
: Statistical test for the existence of a late conservation pattern (high-high-low pattern)
- modules
a list storing three elements for the
ReductiveHourglassTest
,EarlyConservationTest
,LateConservationTest
, orReverseHourglassTest
: early, mid, and late. Each element expects a numeric vector specifying the developmental stages or experiments that correspond to each module. For example:modules
=list(early = 1:2, mid = 3:5, late = 6:7)
divides a dataset storing seven developmental stages into 3 modules.
- permutations
a numeric value specifying the number of permutations to be performed for the
FlatLineTest
,EarlyConservationTest
,LateConservationTest
,ReductiveHourglassTest
orReverseHourglassTest
.- p.value
a boolean value specifying whether the p-value of the test statistic shall be printed within the legend, for each expression set.
- shaded.area
a boolean value specifying whether a shaded area shall be drawn for the developmental stages defined to be the presumptive phylotypic period.
- xlab
label of x-axis.
- ylab
label of y-axis.
- main
figure title.
- legend.title
legend title.
- lwd
line width.
- alpha
transparency of the shaded area and error ribbon (between [0,1]). Default is
alpha = 0.1
.- y.ticks
number of ticks on the y-axis. Default is
ticks = 3
.
Value
a ggplot object visualising the transcriptome index of each given expression set, together with its standard deviation per stage, obtained by permuting the gene ages. The profiles are shown on the same axes, so that they can be readily compared. Optionally, the p-value of each profile, with respect to the choice of statistic, is shown.
Examples
data(PhyloExpressionSetExample)
# remove top 1% expressed genes
genes.top_expression <- TopExpressionGenes(PhyloExpressionSetExample, p=.99)
PhyloExpressionSetExample.top_removed <- subset(PhyloExpressionSetExample,
!(GeneID %in% genes.top_expression))
expression_sets = list(PhyloExpressionSetExample, PhyloExpressionSetExample.top_removed)
set_labels = c("100%", "99%")
# Flat line test
PlotSignatureMultiple(ExpressionSets = expression_sets,
set.labels = set_labels,
TestStatistic="FlatLineTest",
main = "A. thaliana embryogenesis",
legend.title = "Top Expressed Genes Quantile")
#>
#> [ Number of Eigen threads that are employed on your machine: 12 ]
#>
#> [ Computing age assignment permutations for test statistic ... ]
#>
[=========================================] 100%
#> [ Computing variances of permuted transcriptome signatures ... ]
#>
#>
#> Total runtime of your permutation test: 0.128 seconds.
#>
#> -> We recommended using at least 20000 permutations to achieve a sufficient permutation test.
#>
#> [ Number of Eigen threads that are employed on your machine: 12 ]
#>
#> [ Computing age assignment permutations for test statistic ... ]
#>
[=========================================] 100%
#> [ Computing variances of permuted transcriptome signatures ... ]
#>
#>
#> Total runtime of your permutation test: 0.125 seconds.
#>
#> -> We recommended using at least 20000 permutations to achieve a sufficient permutation test.
# Reductive hourglass test
PlotSignatureMultiple(ExpressionSets = expression_sets,
set.labels = set_labels,
TestStatistic="ReductiveHourglassTest",
main = "A. thaliana embryogenesis",
legend.title = "Top Expressed Genes Quantile",
modules=list(early=1:2, mid=3:5, late=6:7),
shaded.area=TRUE)
#>
#> [ Number of Eigen threads that are employed on your machine: 12 ]
#>
#> [ Computing age assignment permutations for test statistic ... ]
#>
[=========================================] 100%
#> [ Computing variances of permuted transcriptome signatures ... ]
#>
#>
#> [ Number of Eigen threads that are employed on your machine: 12 ]
#>
#> [ Computing age assignment permutations for test statistic ... ]
#>
[=========================================] 100%
#> [ Computing variances of permuted transcriptome signatures ... ]
#>