Title: | Tools for data visualization |
---|---|
Description: | systemPipeTools package extends the widely used systemPipeR (SPR) workflow environment with an enhanced toolkit for data visualization, including utilities to automate the data visualizaton for analysis of differentially expressed genes (DEGs). systemPipeTools provides data transformation and data exploration functions via scatterplots, hierarchical clustering heatMaps, principal component analysis, multidimensional scaling, generalized principal components, t-Distributed Stochastic Neighbor embedding (t-SNE), and MA and volcano plots. All these utilities can be integrated with the modular design of the systemPipeR environment that allows users to easily substitute any of these features and/or custom with alternatives. |
Authors: | Daniela Cassol [aut, cre], Ponmathi Ramasamy [aut], Le Zhang [aut], Thomas Girke [aut] |
Maintainer: | Daniela Cassol <[email protected]> |
License: | Artistic-2.0 |
Version: | 1.15.0 |
Built: | 2024-10-31 05:41:21 UTC |
Source: | https://github.com/bioc/systemPipeTools |
Convenience wrapper function to transform raw read counts using
the DESeq2::DESeq2-package()
package transformations methods. The input
file has to contain all the genes, not just differentially expressed ones.
exploreDDS( countMatrix, targets, cmp = cmp[[1]], preFilter = NULL, transformationMethod = "raw", blind = TRUE )
exploreDDS( countMatrix, targets, cmp = cmp[[1]], preFilter = NULL, transformationMethod = "raw", blind = TRUE )
countMatrix |
|
targets |
targets |
cmp |
|
preFilter |
allows removing rows in which there are very few reads.
Accepts a numeric value with the minimum of total reads to keep. Default is
|
transformationMethod |
a |
blind |
logical, whether to blind the transformation to the experimental
design (see varianceStabilizingTransformation),
from |
Note that the recommendation is to use the resulting transformed
values in the transformationMethod
argument only for visualization and
clustering, not for differential expression analysis which needs raw counts.
Users are strongly encouraged to consult the DESeq2::DESeq2-package()
vignette for
more detailed information on this topic and how to properly run DESeq2
on
data sets with more complex experimental designs.
returns an object of class DESeq2::DESeqTransform()
.
Daniela Cassol
For more details on DESeq2
, please consult the following
page:
DESeq2.
For more details on targets
file definition, please consult the following
page:
systemPipeR.
## Targets file targetspath <- system.file("extdata", "targets.txt", package = "systemPipeR") targets <- read.delim(targetspath, comment = "#") cmp <- systemPipeR::readComp(file = targetspath, format = "matrix", delim = "-") ## Count table file countMatrixPath <- system.file("extdata", "countDFeByg.xls", package = "systemPipeR") countMatrix <- read.delim(countMatrixPath, row.names = 1) ## Run exploredds <- exploreDDS(countMatrix, targets, cmp = cmp[[1]], preFilter = NULL, transformationMethod = "raw" ) exploredds
## Targets file targetspath <- system.file("extdata", "targets.txt", package = "systemPipeR") targets <- read.delim(targetspath, comment = "#") cmp <- systemPipeR::readComp(file = targetspath, format = "matrix", delim = "-") ## Count table file countMatrixPath <- system.file("extdata", "countDFeByg.xls", package = "systemPipeR") countMatrix <- read.delim(countMatrixPath, row.names = 1) ## Run exploredds <- exploreDDS(countMatrix, targets, cmp = cmp[[1]], preFilter = NULL, transformationMethod = "raw" ) exploredds
Scatterplot of transformed counts from two samples or grid of all samples
exploreDDSplot( countMatrix, targets, cmp = cmp[[1]], preFilter = NULL, samples, blind = TRUE, scattermatrix = FALSE, plotly = FALSE, savePlot = FALSE, filePlot = NULL )
exploreDDSplot( countMatrix, targets, cmp = cmp[[1]], preFilter = NULL, samples, blind = TRUE, scattermatrix = FALSE, plotly = FALSE, savePlot = FALSE, filePlot = NULL )
countMatrix |
|
targets |
targets |
cmp |
|
preFilter |
allows removing rows in which there are very few reads.
Accepts a numeric value with the minimum of total reads to keep.
Default is |
samples |
a |
blind |
logical, whether to blind the transformation to the experimental
design (see varianceStabilizingTransformation), from |
scattermatrix |
if |
plotly |
logical: when |
savePlot |
logical: when |
filePlot |
file name where the plot will be saved. For more information,
please consult the |
returns an object of ggplot2 plot
.
## Targets file targetspath <- system.file("extdata", "targets.txt", package = "systemPipeR") targets <- read.delim(targetspath, comment = "#") cmp <- systemPipeR::readComp( file = targetspath, format = "matrix", delim = "-" ) ## Count table file countMatrixPath <- system.file("extdata", "countDFeByg.xls", package = "systemPipeR" ) countMatrix <- read.delim(countMatrixPath, row.names = 1) ## Plot exploreDDSplot(countMatrix, targets, cmp = cmp[[1]], preFilter = NULL, samples = c(3, 4) )
## Targets file targetspath <- system.file("extdata", "targets.txt", package = "systemPipeR") targets <- read.delim(targetspath, comment = "#") cmp <- systemPipeR::readComp( file = targetspath, format = "matrix", delim = "-" ) ## Count table file countMatrixPath <- system.file("extdata", "countDFeByg.xls", package = "systemPipeR" ) countMatrix <- read.delim(countMatrixPath, row.names = 1) ## Plot exploreDDSplot(countMatrix, targets, cmp = cmp[[1]], preFilter = NULL, samples = c(3, 4) )
This function computes and plots generalized principal components analysis for dimension reduction of count expression matrix.
GLMplot( exploredds, L = 2, plotly = FALSE, savePlot = FALSE, filePlot = NULL, ... )
GLMplot( exploredds, L = 2, plotly = FALSE, savePlot = FALSE, filePlot = NULL, ... )
exploredds |
object of class |
L |
desired number of latent dimensions (positive integer). |
plotly |
logical: when |
savePlot |
logical: when |
filePlot |
file name where the plot will be saved. For more information,
please consult the |
... |
additional parameters for the |
returns an object of ggplot
or plotly
class.
F. William Townes and Kelly Street (2020). glmpca: Dimension Reduction of Non-Normally Distributed Data. R package version 0.2.0. https://CRAN.R-project.org/package=glmpca
## Targets file targetspath <- system.file("extdata", "targets.txt", package = "systemPipeR") targets <- read.delim(targetspath, comment = "#") cmp <- systemPipeR::readComp(file = targetspath, format = "matrix", delim = "-") ## Count table file countMatrixPath <- system.file("extdata", "countDFeByg.xls", package = "systemPipeR") countMatrix <- read.delim(countMatrixPath, row.names = 1) ## Plot exploredds <- exploreDDS(countMatrix, targets, cmp = cmp[[1]], preFilter = NULL, transformationMethod = "raw") GLMplot(exploredds, plotly = FALSE)
## Targets file targetspath <- system.file("extdata", "targets.txt", package = "systemPipeR") targets <- read.delim(targetspath, comment = "#") cmp <- systemPipeR::readComp(file = targetspath, format = "matrix", delim = "-") ## Count table file countMatrixPath <- system.file("extdata", "countDFeByg.xls", package = "systemPipeR") countMatrix <- read.delim(countMatrixPath, row.names = 1) ## Plot exploredds <- exploreDDS(countMatrix, targets, cmp = cmp[[1]], preFilter = NULL, transformationMethod = "raw") GLMplot(exploredds, plotly = FALSE)
This function computes the sample-wise correlation coefficients
using the stats::cor()
function from the transformed expression values.
After transformation to a distance matrix, hierarchical clustering is
performed with the stats::hclust()
function, and the result is plotted as
a dendrogram.
hclustplot( exploredds, method = "spearman", plotly = FALSE, savePlot = FALSE, filePlot = NULL )
hclustplot( exploredds, method = "spearman", plotly = FALSE, savePlot = FALSE, filePlot = NULL )
exploredds |
object of class |
method |
a |
plotly |
logical: when |
savePlot |
logical: when |
filePlot |
file name where the plot will be saved. For more information,
please consult the |
returns an object of ggplot
or plotly
class.
## Targets file targetspath <- system.file("extdata", "targets.txt", package = "systemPipeR") targets <- read.delim(targetspath, comment = "#") cmp <- systemPipeR::readComp(file = targetspath, format = "matrix", delim = "-") ## Count table file countMatrixPath <- system.file("extdata", "countDFeByg.xls", package = "systemPipeR") countMatrix <- read.delim(countMatrixPath, row.names = 1) ## Plot exploredds <- exploreDDS(countMatrix, targets, cmp = cmp[[1]], preFilter = NULL, transformationMethod = "rlog" ) hclustplot(exploredds, method = "spearman")
## Targets file targetspath <- system.file("extdata", "targets.txt", package = "systemPipeR") targets <- read.delim(targetspath, comment = "#") cmp <- systemPipeR::readComp(file = targetspath, format = "matrix", delim = "-") ## Count table file countMatrixPath <- system.file("extdata", "countDFeByg.xls", package = "systemPipeR") countMatrix <- read.delim(countMatrixPath, row.names = 1) ## Plot exploredds <- exploreDDS(countMatrix, targets, cmp = cmp[[1]], preFilter = NULL, transformationMethod = "rlog" ) hclustplot(exploredds, method = "spearman")
This function performs hierarchical clustering on the transformed expression matrix generated with the DESeq2 package. It uses, by default, a Pearson correlation-based distance measure and complete linkage for cluster join.
heatMaplot( exploredds, clust, DEGlist = NULL, plotly = FALSE, savePlot = FALSE, filePlot = NULL, ... )
heatMaplot( exploredds, clust, DEGlist = NULL, plotly = FALSE, savePlot = FALSE, filePlot = NULL, ... )
exploredds |
object of class |
clust |
select the data to apply the distance matrix computation.
If |
DEGlist |
List of up or down regulated gene/transcript identifiers
meeting the chosen filter settings for all comparisons defined in data
frames |
plotly |
logical: when |
savePlot |
logical: when |
filePlot |
file name where the plot will be saved. For more information,
please consult the |
... |
additional parameters for the |
returns an object of pheatmap
or plotly
class.
Raivo Kolde (2019). pheatmap: Pretty Heatmaps. R package version 1.0.12. https://CRAN.R-project.org/package=pheatmap
### Load data targetspath <- system.file("extdata", "targets.txt", package = "systemPipeR") targets <- read.delim(targetspath, comment = "#") cmp <- systemPipeR::readComp(file = targetspath, format = "matrix", delim = "-") countMatrixPath <- system.file("extdata", "countDFeByg.xls", package = "systemPipeR") countMatrix <- read.delim(countMatrixPath, row.names = 1) ## Samples plot exploredds <- exploreDDS(countMatrix, targets, cmp = cmp[[1]], preFilter = NULL, transformationMethod = "rlog" ) heatMaplot(exploredds, clust = "samples", plotly = TRUE) ## Individuals genes identified in DEG analysis ### DEG analysis with `systemPipeR` degseqDF <- systemPipeR::run_DESeq2( countDF = countMatrix, targets = targets, cmp = cmp[[1]], independent = FALSE ) DEG_list <- systemPipeR::filterDEGs( degDF = degseqDF, filter = c(Fold = 2, FDR = 10) ) ### Plot heatMaplot(exploredds, clust = "ind", DEGlist = unique(as.character(unlist(DEG_list[[1]]))) )
### Load data targetspath <- system.file("extdata", "targets.txt", package = "systemPipeR") targets <- read.delim(targetspath, comment = "#") cmp <- systemPipeR::readComp(file = targetspath, format = "matrix", delim = "-") countMatrixPath <- system.file("extdata", "countDFeByg.xls", package = "systemPipeR") countMatrix <- read.delim(countMatrixPath, row.names = 1) ## Samples plot exploredds <- exploreDDS(countMatrix, targets, cmp = cmp[[1]], preFilter = NULL, transformationMethod = "rlog" ) heatMaplot(exploredds, clust = "samples", plotly = TRUE) ## Individuals genes identified in DEG analysis ### DEG analysis with `systemPipeR` degseqDF <- systemPipeR::run_DESeq2( countDF = countMatrix, targets = targets, cmp = cmp[[1]], independent = FALSE ) DEG_list <- systemPipeR::filterDEGs( degDF = degseqDF, filter = c(Fold = 2, FDR = 10) ) ### Plot heatMaplot(exploredds, clust = "ind", DEGlist = unique(as.character(unlist(DEG_list[[1]]))) )
This function plots log2 fold changes (y-axis) versus the mean of normalized counts (on the x-axis). Statistically significant features are colored.
MAplot( degseqDF, FDR.cutoff = 0.05, comparison, filter = c(Fold = 2, FDR = 10), genes = "NULL", plotly = FALSE, savePlot = FALSE, filePlot = NULL )
MAplot( degseqDF, FDR.cutoff = 0.05, comparison, filter = c(Fold = 2, FDR = 10), genes = "NULL", plotly = FALSE, savePlot = FALSE, filePlot = NULL )
degseqDF |
object of class |
FDR.cutoff |
filter cutoffs for the p-value adjusted. |
comparison |
|
filter |
Named vector with filter cutoffs of format c(Fold=2, FDR=1) where Fold refers to the fold change cutoff (unlogged) and FDR to the p-value cutoff. |
genes |
|
plotly |
logical: when |
savePlot |
logical: when |
filePlot |
file name where the plot will be saved. For more information,
please consult the |
returns an object of ggplot
or plotly
class.
## Load targets file and count reads dataframe targetspath <- system.file("extdata", "targets.txt", package = "systemPipeR") targets <- read.delim(targetspath, comment = "#") cmp <- systemPipeR::readComp( file = targetspath, format = "matrix", delim = "-" ) countMatrixPath <- system.file("extdata", "countDFeByg.xls", package = "systemPipeR" ) countMatrix <- read.delim(countMatrixPath, row.names = 1) ### DEG analysis with `systemPipeR` degseqDF <- systemPipeR::run_DESeq2( countDF = countMatrix, targets = targets, cmp = cmp[[1]], independent = FALSE ) DEG_list <- systemPipeR::filterDEGs( degDF = degseqDF, filter = c(Fold = 2, FDR = 10) ) ## Plot MAplot(degseqDF, comparison = "M12-A12", filter = c(Fold = 1, FDR = 20), genes = "ATCG00280" )
## Load targets file and count reads dataframe targetspath <- system.file("extdata", "targets.txt", package = "systemPipeR") targets <- read.delim(targetspath, comment = "#") cmp <- systemPipeR::readComp( file = targetspath, format = "matrix", delim = "-" ) countMatrixPath <- system.file("extdata", "countDFeByg.xls", package = "systemPipeR" ) countMatrix <- read.delim(countMatrixPath, row.names = 1) ### DEG analysis with `systemPipeR` degseqDF <- systemPipeR::run_DESeq2( countDF = countMatrix, targets = targets, cmp = cmp[[1]], independent = FALSE ) DEG_list <- systemPipeR::filterDEGs( degDF = degseqDF, filter = c(Fold = 2, FDR = 10) ) ## Plot MAplot(degseqDF, comparison = "M12-A12", filter = c(Fold = 1, FDR = 20), genes = "ATCG00280" )
This function computes and plots multidimensional scaling
analysis for dimension reduction of count expression matrix. Internally, it
is applied the stats::dist()
function to the transformed count matrix to
get sample-to-sample distances.
MDSplot( exploredds, method = "spearman", plotly = FALSE, savePlot = FALSE, filePlot = NULL )
MDSplot( exploredds, method = "spearman", plotly = FALSE, savePlot = FALSE, filePlot = NULL )
exploredds |
object of class |
method |
a |
plotly |
logical: when |
savePlot |
logical: when |
filePlot |
file name where the plot will be saved. For more information,
please consult the |
returns an object of ggplot
or plotly
class.
## Targets file targetspath <- system.file("extdata", "targets.txt", package = "systemPipeR") targets <- read.delim(targetspath, comment = "#") cmp <- systemPipeR::readComp(file = targetspath, format = "matrix", delim = "-") ## Count table file countMatrixPath <- system.file("extdata", "countDFeByg.xls", package = "systemPipeR") countMatrix <- read.delim(countMatrixPath, row.names = 1) ## Plot exploredds <- exploreDDS(countMatrix, targets, cmp = cmp[[1]], preFilter = NULL, transformationMethod = "rlog") MDSplot(exploredds, plotly = FALSE)
## Targets file targetspath <- system.file("extdata", "targets.txt", package = "systemPipeR") targets <- read.delim(targetspath, comment = "#") cmp <- systemPipeR::readComp(file = targetspath, format = "matrix", delim = "-") ## Count table file countMatrixPath <- system.file("extdata", "countDFeByg.xls", package = "systemPipeR") countMatrix <- read.delim(countMatrixPath, row.names = 1) ## Plot exploredds <- exploreDDS(countMatrix, targets, cmp = cmp[[1]], preFilter = NULL, transformationMethod = "rlog") MDSplot(exploredds, plotly = FALSE)
This function plots a Principal Component Analysis (PCA) from transformed expression matrix. This plot shows samples variation based on the expression values and identifies batch effects.
PCAplot(exploredds, plotly = FALSE, savePlot = FALSE, filePlot = NULL)
PCAplot(exploredds, plotly = FALSE, savePlot = FALSE, filePlot = NULL)
exploredds |
object of class |
plotly |
logical: when |
savePlot |
logical: when |
filePlot |
file name where the plot will be saved. For more information,
please consult the |
returns an object of ggplot
or plotly
class.
## Targets file targetspath <- system.file("extdata", "targets.txt", package = "systemPipeR") targets <- read.delim(targetspath, comment = "#") cmp <- systemPipeR::readComp(file = targetspath, format = "matrix", delim = "-") ## Count table file countMatrixPath <- system.file("extdata", "countDFeByg.xls", package = "systemPipeR") countMatrix <- read.delim(countMatrixPath, row.names = 1) ## Plot exploredds <- exploreDDS(countMatrix, targets, cmp = cmp[[1]], preFilter = NULL, transformationMethod = "rlog") PCAplot(exploredds, plotly = TRUE)
## Targets file targetspath <- system.file("extdata", "targets.txt", package = "systemPipeR") targets <- read.delim(targetspath, comment = "#") cmp <- systemPipeR::readComp(file = targetspath, format = "matrix", delim = "-") ## Count table file countMatrixPath <- system.file("extdata", "countDFeByg.xls", package = "systemPipeR") countMatrix <- read.delim(countMatrixPath, row.names = 1) ## Plot exploredds <- exploreDDS(countMatrix, targets, cmp = cmp[[1]], preFilter = NULL, transformationMethod = "rlog") PCAplot(exploredds, plotly = TRUE)
Create an HTML table using DT package with fixed columns
showDT(data, ...)
showDT(data, ...)
data |
data object (either a matrix or a data frame). |
... |
Additional arguments used by dDT::atatable() function. |
returns an object of datatables
and htmlwidget
.
showDT(iris)
showDT(iris)
This function computes and plots t-Distributed Stochastic
Neighbor embedding (t-SNE) analysis for unsupervised nonlinear
dimensionality reduction of count expression matrix. Internally, it is
applied the
Rtsne::Rtsne()
function, using the exact t-SNE computing with theta=0.0
.
tSNEplot( countMatrix, targets, plotly = FALSE, savePlot = FALSE, filePlot = NULL, ... )
tSNEplot( countMatrix, targets, plotly = FALSE, savePlot = FALSE, filePlot = NULL, ... )
countMatrix |
|
targets |
targets |
plotly |
logical: when |
savePlot |
logical: when |
filePlot |
file name where the plot will be saved. For more
information, please consult the |
... |
additional parameters for the |
returns an object of ggplot
or plotly
class.
Jesse H. Krijthe (2015). Rtsne: T-Distributed Stochastic Neighbor Embedding using a Barnes-Hut Implementation, URL: https://github.com/jkrijthe/Rtsne
targetspath <- system.file("extdata", "targets.txt", package = "systemPipeR") targets <- read.delim(targetspath, comment = "#") cmp <- systemPipeR::readComp(file = targetspath, format = "matrix", delim = "-") countMatrixPath <- system.file("extdata", "countDFeByg.xls", package = "systemPipeR") countMatrix <- read.delim(countMatrixPath, row.names = 1) set.seed(42) tSNEplot(countMatrix, targets, perplexity = 5)
targetspath <- system.file("extdata", "targets.txt", package = "systemPipeR") targets <- read.delim(targetspath, comment = "#") cmp <- systemPipeR::readComp(file = targetspath, format = "matrix", delim = "-") countMatrixPath <- system.file("extdata", "countDFeByg.xls", package = "systemPipeR") countMatrix <- read.delim(countMatrixPath, row.names = 1) set.seed(42) tSNEplot(countMatrix, targets, perplexity = 5)
volcanoplot
A simple function that shows statistical significance
(p-value
) versus magnitude of change (log2 fold change
).
volcanoplot( degseqDF, comparison, filter = c(Fold = 2, FDR = 10), genes = "NULL", plotly = FALSE, savePlot = FALSE, filePlot = NULL )
volcanoplot( degseqDF, comparison, filter = c(Fold = 2, FDR = 10), genes = "NULL", plotly = FALSE, savePlot = FALSE, filePlot = NULL )
degseqDF |
object of class |
comparison |
|
filter |
Named vector with filter cutoffs of format c(Fold=2, FDR=1) where Fold refers to the fold change cutoff (unlogged) and FDR to the p-value cutoff. |
genes |
|
plotly |
logical: when |
savePlot |
logical: when |
filePlot |
file name where the plot will be saved. For more information,
please consult the |
returns an object of ggplot
or plotly
class.
## Load targets file and count reads dataframe targetspath <- system.file("extdata", "targets.txt", package = "systemPipeR") targets <- read.delim(targetspath, comment = "#") cmp <- systemPipeR::readComp( file = targetspath, format = "matrix", delim = "-" ) countMatrixPath <- system.file("extdata", "countDFeByg.xls", package = "systemPipeR" ) countMatrix <- read.delim(countMatrixPath, row.names = 1) ### DEG analysis with `systemPipeR` degseqDF <- systemPipeR::run_DESeq2( countDF = countMatrix, targets = targets, cmp = cmp[[1]], independent = FALSE) DEG_list <- systemPipeR::filterDEGs( degDF = degseqDF, filter = c(Fold = 2, FDR = 10)) ## Plot volcanoplot(degseqDF, comparison = "M12-A12", filter = c(Fold = 1, FDR = 20), genes = "ATCG00280")
## Load targets file and count reads dataframe targetspath <- system.file("extdata", "targets.txt", package = "systemPipeR") targets <- read.delim(targetspath, comment = "#") cmp <- systemPipeR::readComp( file = targetspath, format = "matrix", delim = "-" ) countMatrixPath <- system.file("extdata", "countDFeByg.xls", package = "systemPipeR" ) countMatrix <- read.delim(countMatrixPath, row.names = 1) ### DEG analysis with `systemPipeR` degseqDF <- systemPipeR::run_DESeq2( countDF = countMatrix, targets = targets, cmp = cmp[[1]], independent = FALSE) DEG_list <- systemPipeR::filterDEGs( degDF = degseqDF, filter = c(Fold = 2, FDR = 10)) ## Plot volcanoplot(degseqDF, comparison = "M12-A12", filter = c(Fold = 1, FDR = 20), genes = "ATCG00280")