1 Getting started

Rvisdiff is an R package distributed as part of the Bioconductor project. To install the package, start R and enter:

if (!requireNamespace("BiocManager", quietly = TRUE)) {
    install.packages("BiocManager")
}
BiocManager::install("Rvisdiff")

The GitHub repository for Rvisdiff is https://github.com/BioinfoUSAL/Rvisdiff. This is the place to file an issue, report a bug, or provide a pull request.

Once Rvisdiff is installed, it can be loaded by the following command.

library("Rvisdiff")

2 Introduction

Differential expression analysis generates a big report which needs a manual inspection for the optimization and interpretation of results. Researchers have designed visualization techniques to facilitate these tasks but their generation with code or statistics packages avoids the quick and massive exploration of results. We have designed Rvisdiff to integrate graphs in an easy to use and interactive web page.The user can explore the differential expression results and the source expression data in the same view.

As input data the package receives two tables with the differential expression results and the raw/normalized expression data. It detects the default output of DESeq2, edgeR and limma packages and no data conversion is needed. The user can also generate a custom data frame which integrates a statistical testing output with a fold change and mean calculation for each variable.

As output the package generates a local HTML page that can be seen in a Web browser. It is not necessary the installation of additional software such as application servers or programming languages. This feature ensures portability and ease of use. Moreover, results are stored in the local computer, avoiding any network sharing or data upload to external servers, which ensures the data privacy.

3 Input data

In this example we use as input the airway data package which contains the read counts in genes for an RNA-Seq experiment on four human airway smooth muscle cell lines treated with dexamethasone. The code below shows how to load the package and the data extraction of main data features that we need for the differential expression analysis and the posterior visualization with Rvisdiff. The countdata variable contains a data frame with the number of sequence counts for each gene (rows) and sample (columns). The coldata variable contains input phenotypes for the differential expression analysis and its posterior representation.

The following code loads the necessary libraries and formats the input sample conditions.

library(Rvisdiff)
library(airway)
data("airway")
se <- airway
se$dex <- relevel(se$dex, ref="untrt")
countdata <- assay(se)
coldata <- colData(se)

4 Generating the Report

4.1 Generating Report From DESeq2 results

The code below shows how to perform a differential expression analysis with DESeq2 and its representation with Rvisdiff.

library(DESeq2)
dds <- DESeqDataSet(se, design = ~ cell + dex)
dds <- DESeq(dds)
dres <- results(dds, independentFiltering = FALSE)
DEreport(dres, countdata, coldata$dex)

4.2 Generating Report From edgeR results

The code below shows how to perform a differential expression analysis with edgeR and its representation with Rvisdiff.

library(edgeR)
design <- model.matrix(~ cell + dex, data = coldata)
dl <- DGEList(counts = countdata, group = coldata$dex)
dl <- calcNormFactors(dl)
dl <- estimateDisp(dl, design=design)
de <- exactTest(dl,pair=1:2)
tt <- topTags(de, n = Inf, adjust.method = "BH", sort.by = "none")
DEreport(tt, countdata, coldata$dex) 

4.3 Generating Report From limma results

The code below shows how to perform a differential expression analysis with limma and its representation with Rvisdiff.

library(limma)
design <- model.matrix(~ 0 + dex + cell, data = coldata)
contr <- makeContrasts(dextrt - dexuntrt,levels=colnames(design))
limmaexprs <- voom(countdata, design)
fit <- lmFit(limmaexprs, design)
fit <- contrasts.fit(fit, contrasts=contr)
fit <- eBayes(fit)
limmares <- topTable(fit, coef = 1, number = Inf, sort.by = "none",
    adjust.method = "BH")
DEreport(limmares, countdata, coldata$dex) 

4.4 Generating Report From Differential test results

The code below shows how to perform a Wilcoxon test with expression data and its representation with Rvisdiff. This example can be also followed for the representation of resulting analysis from differential means tests.

untrt <- countdata[,coldata$dex=="untrt"]
trt <- countdata[,coldata$dex=="trt"]

library(matrixTests)
wilcox <- col_wilcoxon_twosample(t(untrt), t(trt))
stat <- wilcox$statistic
p <- wilcox$pvalue
log2FoldChange <- log2(rowMeans(trt)+1) - log2(rowMeans(untrt)+1)
wilcox <- cbind(stat = stat, pvalue = round(p, 6),
    padj = p.adjust(wilcox[,2], method = "BH"),
    baseMean = rowMeans(countdata),
    log2FoldChange = log2FoldChange)
rownames(wilcox) <- rownames(countdata)

DEreport(wilcox, countdata, coldata$dex)

5 Resulting Graphical User Interface

Figure 1 shows the resulting Web page generated with the DEreport function. The user can select which genes appear in the graphs selecting them in the results table. It contains the following graphs:

  • Volcano Plot: It is a scatter plot in which the values of rate of change are plotted in logarithmic scale (log2foldchange) versus the p-value resulting from the contrast test is scale minus logarithm 10 (-log10pvalue). Points are highlighted when the mouse is hovered over the results table. Variable name appears on point mouse over.
  • MA-Plot: is a scatter plot showing mean expression values versus rate of change, both are plotted in logarithmic scale to avoid excessive scatter. It has the same interactivity features as Volcano plots.
  • Line diagram: the gene expression levels (ordinates) in each sample (abscissae) are represented as a line. Diagram is divided based on input phenotype.
  • Box plot: they allow us to visualize the distribution, degree of asymmetry, extreme values and value of the median. It is also useful for comparing two distributions if we represent them in the same graph. The resulting graphs show the difference in expression between genes or conditions.
  • Cluster Heatmap: expression data are displayed in a grid where each row represents a gene and each column represents a sample. The color and intensity of the boxes are used to represent changes (usually scaled per gene, avoiding absolute values) in gene expression. The heatmap shows also a clustering tree that groups genes and samples based on the similarity of their gene expression pattern. The user can change the color scale and toggle rendering from raw to scaled values. Moreover, the graph provides a zoom feature which enables to set the focus on a set of samples or genes.

Figure 1 web interface

SessionInfo

sessionInfo()
## R version 4.4.0 beta (2024-04-15 r86425)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 22.04.4 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.19-bioc/R/lib/libRblas.so 
## LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.10.0
## 
## locale:
##  [1] LC_CTYPE=en_US.UTF-8       LC_NUMERIC=C              
##  [3] LC_TIME=en_US.UTF-8        LC_COLLATE=en_US.UTF-8    
##  [5] LC_MONETARY=en_US.UTF-8    LC_MESSAGES=en_US.UTF-8   
##  [7] LC_PAPER=en_US.UTF-8       LC_NAME=C                 
##  [9] LC_ADDRESS=C               LC_TELEPHONE=C            
## [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C       
## 
## time zone: America/New_York
## tzcode source: system (glibc)
## 
## attached base packages:
## [1] stats4    stats     graphics  grDevices utils     datasets  methods  
## [8] base     
## 
## other attached packages:
##  [1] matrixTests_0.2.3           edgeR_4.2.0                
##  [3] limma_3.60.0                DESeq2_1.44.0              
##  [5] airway_1.23.0               SummarizedExperiment_1.34.0
##  [7] Biobase_2.64.0              GenomicRanges_1.56.0       
##  [9] GenomeInfoDb_1.40.0         IRanges_2.38.0             
## [11] S4Vectors_0.42.0            BiocGenerics_0.50.0        
## [13] MatrixGenerics_1.16.0       matrixStats_1.3.0          
## [15] Rvisdiff_1.2.0              BiocStyle_2.32.0           
## 
## loaded via a namespace (and not attached):
##  [1] gtable_0.3.5            xfun_0.43               bslib_0.7.0            
##  [4] ggplot2_3.5.1           lattice_0.22-6          vctrs_0.6.5            
##  [7] tools_4.4.0             generics_0.1.3          parallel_4.4.0         
## [10] tibble_3.2.1            fansi_1.0.6             highr_0.10             
## [13] pkgconfig_2.0.3         Matrix_1.7-0            lifecycle_1.0.4        
## [16] GenomeInfoDbData_1.2.12 compiler_4.4.0          statmod_1.5.0          
## [19] munsell_0.5.1           codetools_0.2-20        htmltools_0.5.8.1      
## [22] sass_0.4.9              yaml_2.3.8              pillar_1.9.0           
## [25] crayon_1.5.2            jquerylib_0.1.4         BiocParallel_1.38.0    
## [28] DelayedArray_0.30.0     cachem_1.0.8            abind_1.4-5            
## [31] tidyselect_1.2.1        locfit_1.5-9.9          digest_0.6.35          
## [34] dplyr_1.1.4             bookdown_0.39           splines_4.4.0          
## [37] fastmap_1.1.1           grid_4.4.0              colorspace_2.1-0       
## [40] cli_3.6.2               SparseArray_1.4.0       magrittr_2.0.3         
## [43] S4Arrays_1.4.0          utf8_1.2.4              scales_1.3.0           
## [46] UCSC.utils_1.0.0        rmarkdown_2.26          XVector_0.44.0         
## [49] httr_1.4.7              evaluate_0.23           knitr_1.46             
## [52] rlang_1.1.3             Rcpp_1.0.12             glue_1.7.0             
## [55] BiocManager_1.30.22     jsonlite_1.8.8          R6_2.5.1               
## [58] zlibbioc_1.50.0