Rvisdiff is an R/Bioconductor package which generates an interactive interface for the interpretation of differential expression results. It generates a local Web page which enables the exploration of statistical analysis results with the generation of auto-analytical visualizations. The package supports as input the output of popular differential expression packages such as DESeq2, edgeR and limma.
Rvisdiff 1.6.0
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")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.
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)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)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) 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) 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)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:
Figure 1 web interface
sessionInfo()## R version 4.5.0 RC (2025-04-04 r88126)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.2 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.21-bioc/R/lib/libRblas.so 
## LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.12.0  LAPACK version 3.12.0
## 
## locale:
##  [1] LC_CTYPE=en_US.UTF-8       LC_NUMERIC=C              
##  [3] LC_TIME=en_GB              LC_COLLATE=C              
##  [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.6.0                
##  [3] limma_3.64.0                DESeq2_1.48.0              
##  [5] airway_1.27.0               SummarizedExperiment_1.38.0
##  [7] Biobase_2.68.0              GenomicRanges_1.60.0       
##  [9] GenomeInfoDb_1.44.0         IRanges_2.42.0             
## [11] S4Vectors_0.46.0            BiocGenerics_0.54.0        
## [13] generics_0.1.3              MatrixGenerics_1.20.0      
## [15] matrixStats_1.5.0           Rvisdiff_1.6.0             
## [17] BiocStyle_2.36.0           
## 
## loaded via a namespace (and not attached):
##  [1] gtable_0.3.6            xfun_0.52               bslib_0.9.0            
##  [4] ggplot2_3.5.2           lattice_0.22-7          vctrs_0.6.5            
##  [7] tools_4.5.0             parallel_4.5.0          tibble_3.2.1           
## [10] pkgconfig_2.0.3         Matrix_1.7-3            lifecycle_1.0.4        
## [13] GenomeInfoDbData_1.2.14 compiler_4.5.0          statmod_1.5.0          
## [16] munsell_0.5.1           codetools_0.2-20        htmltools_0.5.8.1      
## [19] sass_0.4.10             yaml_2.3.10             pillar_1.10.2          
## [22] crayon_1.5.3            jquerylib_0.1.4         BiocParallel_1.42.0    
## [25] DelayedArray_0.34.0     cachem_1.1.0            abind_1.4-8            
## [28] tidyselect_1.2.1        locfit_1.5-9.12         digest_0.6.37          
## [31] dplyr_1.1.4             bookdown_0.43           fastmap_1.2.0          
## [34] grid_4.5.0              colorspace_2.1-1        cli_3.6.4              
## [37] SparseArray_1.8.0       magrittr_2.0.3          S4Arrays_1.8.0         
## [40] scales_1.3.0            UCSC.utils_1.4.0        rmarkdown_2.29         
## [43] XVector_0.48.0          httr_1.4.7              evaluate_1.0.3         
## [46] knitr_1.50              rlang_1.1.6             Rcpp_1.0.14            
## [49] glue_1.8.0              BiocManager_1.30.25     jsonlite_2.0.0         
## [52] R6_2.6.1