as.<class>(object), where <class> is the name of the class to which object should be coerced. Additionally, DEFormats extends the original packages by offering constructors for common count data input formats, e.g. SummarizedExperiment objects.
DEFormats 1.20.0
The following examples demonstrate how to shuttle data between DESeqDataSet and DGEList containers. Before we start, let’s first load the required packages.
library(DESeq2)
library(edgeR)
library(DEFormats)We will illustrate the functionality of the package on a mock expression data set of an RNA-seq experiment. The sample counts table can be generated using the function provided by DEFormats:
counts = simulateRnaSeqData()The resulting object is a matrix with rows corresponding to genomic features and columns to samples.
head(counts)##       sample1 sample2 sample3 sample4 sample5 sample6
## gene1      85      76     103     107     140     124
## gene2       1       6      11       6       1       4
## gene3      80      98      39      82      97     113
## gene4      92      83      59      85     100      98
## gene5      36      24      18      50      22      15
## gene6       0       0       1       4       2       3In order to construct a DGEList object we need to provide, apart from the counts matrix, the sample grouping.
group = rep(c("A", "B"), each = 3)
dge = DGEList(counts, group = group)
dge## An object of class "DGEList"
## $counts
##       sample1 sample2 sample3 sample4 sample5 sample6
## gene1      85      76     103     107     140     124
## gene2       1       6      11       6       1       4
## gene3      80      98      39      82      97     113
## gene4      92      83      59      85     100      98
## gene5      36      24      18      50      22      15
## 995 more rows ...
## 
## $samples
##         group lib.size norm.factors
## sample1     A    42832            1
## sample2     A    40511            1
## sample3     A    39299            1
## sample4     B    43451            1
## sample5     B    40613            1
## sample6     B    43662            1A DGEList object can be easily converted to a DESeqDataSet object with the
help of the function as.DESeqDataSet.
dds = as.DESeqDataSet(dge)
dds## class: DESeqDataSet 
## dim: 1000 6 
## metadata(1): version
## assays(1): counts
## rownames(1000): gene1 gene2 ... gene999 gene1000
## rowData names(0):
## colnames(6): sample1 sample2 ... sample5 sample6
## colData names(3): group lib.size norm.factorsJust to make sure that the coercions conserve the data and metadata, we now
convert dds back to a DGEList and compare the result to the original dge
object.
identical(dge, as.DGEList(dds))## [1] TRUENote that because of the use of reference classes in the SummarizedExperiment
class which DESeqDataSet extends, identical will return FALSE for any two
DESeqDataSet instances, even for ones constructed from the same input:
dds1 = DESeqDataSetFromMatrix(counts, data.frame(condition=group), ~ condition)## Warning in DESeqDataSet(se, design = design, ignoreRank): some variables in
## design formula are characters, converting to factorsdds2 = DESeqDataSetFromMatrix(counts, data.frame(condition=group), ~ condition)## Warning in DESeqDataSet(se, design = design, ignoreRank): some variables in
## design formula are characters, converting to factorsidentical(dds1, dds2)## [1] TRUEInstead of a count matrix, simulateRnaSeqData can also return an
annotated RangedSummarizedExperiment object. The advantage of such an object
is that, apart from the counts matrix stored in the assay slot, it also
contains sample description in colData, and gene information stored in
rowRanges as a GRanges object.
se = simulateRnaSeqData(output = "RangedSummarizedExperiment")
se## class: RangedSummarizedExperiment 
## dim: 1000 6 
## metadata(1): version
## assays(1): counts
## rownames(1000): gene1 gene2 ... gene999 gene1000
## rowData names(3): trueIntercept trueBeta trueDisp
## colnames(6): sample1 sample2 ... sample5 sample6
## colData names(1): conditionThe se object can be readily used to construct a DESeqDataSet object,
dds = DESeqDataSet(se, design = ~ condition)which can be converted to a DGEList using the familiar method.
dge = as.DGEList(dds)
dge## An object of class "DGEList"
## $counts
##       sample1 sample2 sample3 sample4 sample5 sample6
## gene1      85      76     103     107     140     124
## gene2       1       6      11       6       1       4
## gene3      80      98      39      82      97     113
## gene4      92      83      59      85     100      98
## gene5      36      24      18      50      22      15
## 995 more rows ...
## 
## $samples
##         group lib.size norm.factors
## sample1     A    42832            1
## sample2     A    40511            1
## sample3     A    39299            1
## sample4     B    43451            1
## sample5     B    40613            1
## sample6     B    43662            1
## 
## $genes
##       seqnames start end width strand trueIntercept trueBeta  trueDisp
## gene1        1     1 100   100      *      6.525909        0 0.1434076
## gene2        1   101 200   100      *      3.347533        0 0.4929634
## gene3        1   201 300   100      *      6.659599        0 0.1395659
## gene4        1   301 400   100      *      6.544859        0 0.1428412
## gene5        1   401 500   100      *      4.829283        0 0.2407022
## 995 more rows ...Note the additional genes element on the dge list compared to the object
from the previous section.
Similarly as for the DESeqDataSet constructor from DESeq2,
it is possible to directly use RangedSummarizedExperiment objects as input to
the generic DGEList constructor defined in DEFormats. This
allows you to use common input objects regardless of whether you are applying
DESeq2 or edgeR to perform your analysis, or to
easily switch between these two tools in your pipeline.
names(colData(se)) = "group"
dge = DGEList(se)
dge## An object of class "DGEList"
## $counts
##       sample1 sample2 sample3 sample4 sample5 sample6
## gene1      85      76     103     107     140     124
## gene2       1       6      11       6       1       4
## gene3      80      98      39      82      97     113
## gene4      92      83      59      85     100      98
## gene5      36      24      18      50      22      15
## 995 more rows ...
## 
## $samples
##         group lib.size norm.factors
## sample1     A    42832            1
## sample2     A    40511            1
## sample3     A    39299            1
## sample4     B    43451            1
## sample5     B    40613            1
## sample6     B    43662            1
## 
## $genes
##       seqnames start end width strand trueIntercept trueBeta  trueDisp
## gene1        1     1 100   100      *      6.525909        0 0.1434076
## gene2        1   101 200   100      *      3.347533        0 0.4929634
## gene3        1   201 300   100      *      6.659599        0 0.1395659
## gene4        1   301 400   100      *      6.544859        0 0.1428412
## gene5        1   401 500   100      *      4.829283        0 0.2407022
## 995 more rows ...We renamed the condition column of colData(se) to group in order to allow
edgeR to automatically pick up the correct sample grouping.
Another way of specifying this is through the argument group to DGEList.
Even though there is no direct method to go from a DGEList to a RangedSummarizedExperiment, the coercion can be easily performed by first converting the object to a DESeqDataSet, which is a subclass of RangedSummarizedExperiment, and then coercing the resulting object to its superclass.
dds = as.DESeqDataSet(dge)
rse = as(dds, "RangedSummarizedExperiment")
rse## class: RangedSummarizedExperiment 
## dim: 1000 6 
## metadata(1): version
## assays(1): counts
## rownames(1000): gene1 gene2 ... gene999 gene1000
## rowData names(3): trueIntercept trueBeta trueDisp
## colnames(6): sample1 sample2 ... sample5 sample6
## colData names(3): group lib.size norm.factorssessionInfo()## R version 4.1.0 (2021-05-18)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 20.04.2 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.13-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.13-bioc/R/lib/libRlapack.so
## 
## 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       
## 
## attached base packages:
## [1] parallel  stats4    stats     graphics  grDevices utils     datasets 
## [8] methods   base     
## 
## other attached packages:
##  [1] DEFormats_1.20.0            edgeR_3.34.0               
##  [3] limma_3.48.0                DESeq2_1.32.0              
##  [5] SummarizedExperiment_1.22.0 Biobase_2.52.0             
##  [7] MatrixGenerics_1.4.0        matrixStats_0.58.0         
##  [9] GenomicRanges_1.44.0        GenomeInfoDb_1.28.0        
## [11] IRanges_2.26.0              S4Vectors_0.30.0           
## [13] BiocGenerics_0.38.0         BiocStyle_2.20.0           
## 
## loaded via a namespace (and not attached):
##  [1] bitops_1.0-7           bit64_4.0.5            RColorBrewer_1.1-2    
##  [4] httr_1.4.2             tools_4.1.0            backports_1.2.1       
##  [7] bslib_0.2.5.1          utf8_1.2.1             R6_2.5.0              
## [10] DBI_1.1.1              colorspace_2.0-1       tidyselect_1.1.1      
## [13] bit_4.0.4              compiler_4.1.0         DelayedArray_0.18.0   
## [16] bookdown_0.22          sass_0.4.0             scales_1.1.1          
## [19] checkmate_2.0.0        genefilter_1.74.0      stringr_1.4.0         
## [22] digest_0.6.27          rmarkdown_2.8          XVector_0.32.0        
## [25] pkgconfig_2.0.3        htmltools_0.5.1.1      fastmap_1.1.0         
## [28] rlang_0.4.11           RSQLite_2.2.7          jquerylib_0.1.4       
## [31] generics_0.1.0         jsonlite_1.7.2         BiocParallel_1.26.0   
## [34] dplyr_1.0.6            RCurl_1.98-1.3         magrittr_2.0.1        
## [37] GenomeInfoDbData_1.2.6 Matrix_1.3-3           Rcpp_1.0.6            
## [40] munsell_0.5.0          fansi_0.4.2            lifecycle_1.0.0       
## [43] stringi_1.6.2          yaml_2.2.1             zlibbioc_1.38.0       
## [46] grid_4.1.0             blob_1.2.1             crayon_1.4.1          
## [49] lattice_0.20-44        Biostrings_2.60.0      splines_4.1.0         
## [52] annotate_1.70.0        KEGGREST_1.32.0        locfit_1.5-9.4        
## [55] knitr_1.33             pillar_1.6.1           geneplotter_1.70.0    
## [58] XML_3.99-0.6           glue_1.4.2             evaluate_0.14         
## [61] data.table_1.14.0      BiocManager_1.30.15    png_0.1-7             
## [64] vctrs_0.3.8            gtable_0.3.0           purrr_0.3.4           
## [67] assertthat_0.2.1       cachem_1.0.5           ggplot2_3.3.3         
## [70] xfun_0.23              xtable_1.8-4           survival_3.2-11       
## [73] tibble_3.1.2           AnnotationDbi_1.54.0   memoise_2.0.0         
## [76] ellipsis_0.3.2