Contents

1 Overview

The alabaster.ranges package implements methods to save genomic ranges (i.e., GRanges and GRangesList objects) to file artifacts and load them back into R. It also supports various CompressedList subclasses, including the somewhat useful CompressedSplitDataFrameList. Check out alabaster.base for more details on the motivation and concepts of the alabaster framework.

2 Quick start

Given some genomic ranges, we can use saveObject() to save it inside a staging directory:

library(GenomicRanges)
gr <- GRanges("chrA", IRanges(sample(100), width=sample(100)))
mcols(gr)$score <- runif(length(gr))
metadata(gr)$genome <- "Aaron"
seqlengths(gr) <- c(chrA=1000)

library(alabaster.ranges)
tmp <- tempfile()
saveObject(gr, tmp)

list.files(tmp, recursive=TRUE)
## [1] "OBJECT"                                 
## [2] "other_annotations/OBJECT"               
## [3] "other_annotations/list_contents.json.gz"
## [4] "range_annotations/OBJECT"               
## [5] "range_annotations/basic_columns.h5"     
## [6] "ranges.h5"                              
## [7] "sequence_information/OBJECT"            
## [8] "sequence_information/info.h5"

We can then easily load it back in with readObject().

roundtrip <- readObject(tmp)
roundtrip
## GRanges object with 100 ranges and 1 metadata column:
##         seqnames    ranges strand |     score
##            <Rle> <IRanges>  <Rle> | <numeric>
##     [1]     chrA    36-121      * |  0.679514
##     [2]     chrA     41-65      * |  0.120720
##     [3]     chrA      2-92      * |  0.555709
##     [4]     chrA     31-64      * |  0.875769
##     [5]     chrA    86-123      * |  0.732745
##     ...      ...       ...    ... .       ...
##    [96]     chrA     29-44      * |  0.970290
##    [97]     chrA    91-167      * |  0.373139
##    [98]     chrA     38-70      * |  0.645407
##    [99]     chrA     14-62      * |  0.463452
##   [100]     chrA    22-113      * |  0.692163
##   -------
##   seqinfo: 1 sequence from an unspecified genome

The same can be done for GRangesList and CompressedList subclasses.

2.1 Further comments

Metadata is preserved during this round-trip:

metadata(roundtrip)
## $genome
## [1] "Aaron"
mcols(roundtrip)
## DataFrame with 100 rows and 1 column
##         score
##     <numeric>
## 1    0.679514
## 2    0.120720
## 3    0.555709
## 4    0.875769
## 5    0.732745
## ...       ...
## 96   0.970290
## 97   0.373139
## 98   0.645407
## 99   0.463452
## 100  0.692163
seqinfo(roundtrip)
## Seqinfo object with 1 sequence from an unspecified genome:
##   seqnames seqlengths isCircular genome
##   chrA           1000         NA   <NA>

Session information

sessionInfo()
## R Under development (unstable) (2024-11-20 r87352)
## Platform: x86_64-apple-darwin20
## Running under: macOS Monterey 12.7.6
## 
## Matrix products: default
## BLAS:   /Library/Frameworks/R.framework/Versions/4.5-x86_64/Resources/lib/libRblas.0.dylib 
## LAPACK: /Library/Frameworks/R.framework/Versions/4.5-x86_64/Resources/lib/libRlapack.dylib;  LAPACK version 3.12.0
## 
## locale:
## [1] C/en_US.UTF-8/en_US.UTF-8/C/en_US.UTF-8/en_US.UTF-8
## 
## time zone: America/New_York
## tzcode source: internal
## 
## attached base packages:
## [1] stats4    stats     graphics  grDevices utils     datasets  methods  
## [8] base     
## 
## other attached packages:
## [1] alabaster.ranges_1.7.0 alabaster.base_1.7.2   GenomicRanges_1.59.1  
## [4] GenomeInfoDb_1.43.1    IRanges_2.41.1         S4Vectors_0.45.2      
## [7] BiocGenerics_0.53.3    generics_0.1.3         BiocStyle_2.35.0      
## 
## loaded via a namespace (and not attached):
##  [1] httr_1.4.7              cli_3.6.3               knitr_1.49             
##  [4] rlang_1.1.4             xfun_0.49               UCSC.utils_1.3.0       
##  [7] jsonlite_1.8.9          htmltools_0.5.8.1       sass_0.4.9             
## [10] rmarkdown_2.29          evaluate_1.0.1          jquerylib_0.1.4        
## [13] fastmap_1.2.0           Rhdf5lib_1.29.0         alabaster.schemas_1.7.0
## [16] yaml_2.3.10             lifecycle_1.0.4         bookdown_0.41          
## [19] BiocManager_1.30.25     compiler_4.5.0          Rcpp_1.0.13-1          
## [22] rhdf5filters_1.19.0     XVector_0.47.0          rhdf5_2.51.0           
## [25] digest_0.6.37           R6_2.5.1                GenomeInfoDbData_1.2.13
## [28] bslib_0.8.0             tools_4.5.0             zlibbioc_1.53.0        
## [31] cachem_1.1.0