alabaster.ranges 1.9.0
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.
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] "_environment.json"
## [3] "other_annotations/OBJECT"
## [4] "other_annotations/list_contents.json.gz"
## [5] "range_annotations/OBJECT"
## [6] "range_annotations/basic_columns.h5"
## [7] "ranges.h5"
## [8] "sequence_information/OBJECT"
## [9] "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 93-162 * | 0.460814
## [2] chrA 32-92 * | 0.988061
## [3] chrA 29-54 * | 0.548615
## [4] chrA 53-79 * | 0.396294
## [5] chrA 56-59 * | 0.254640
## ... ... ... ... . ...
## [96] chrA 96-119 * | 0.901308
## [97] chrA 72-146 * | 0.626908
## [98] chrA 36-131 * | 0.635013
## [99] chrA 91-181 * | 0.646297
## [100] chrA 57-72 * | 0.630350
## -------
## seqinfo: 1 sequence from an unspecified genome
The same can be done for GRangesList
and CompressedList
subclasses.
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.460814
## 2 0.988061
## 3 0.548615
## 4 0.396294
## 5 0.254640
## ... ...
## 96 0.901308
## 97 0.626908
## 98 0.635013
## 99 0.646297
## 100 0.630350
seqinfo(roundtrip)
## Seqinfo object with 1 sequence from an unspecified genome:
## seqnames seqlengths isCircular genome
## chrA 1000 NA <NA>
sessionInfo()
## R version 4.5.0 beta (2025-04-02 r88102)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.2 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.22-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] alabaster.ranges_1.9.0 alabaster.base_1.9.0 GenomicRanges_1.61.0
## [4] GenomeInfoDb_1.45.0 IRanges_2.43.0 S4Vectors_0.47.0
## [7] BiocGenerics_0.55.0 generics_0.1.3 BiocStyle_2.37.0
##
## loaded via a namespace (and not attached):
## [1] httr_1.4.7 cli_3.6.4 knitr_1.50
## [4] rlang_1.1.6 xfun_0.52 UCSC.utils_1.5.0
## [7] jsonlite_2.0.0 htmltools_0.5.8.1 sass_0.4.10
## [10] rmarkdown_2.29 evaluate_1.0.3 jquerylib_0.1.4
## [13] fastmap_1.2.0 Rhdf5lib_1.31.0 alabaster.schemas_1.9.0
## [16] yaml_2.3.10 lifecycle_1.0.4 bookdown_0.43
## [19] BiocManager_1.30.25 compiler_4.5.0 Rcpp_1.0.14
## [22] rhdf5filters_1.21.0 XVector_0.49.0 rhdf5_2.53.0
## [25] digest_0.6.37 R6_2.6.1 GenomeInfoDbData_1.2.14
## [28] bslib_0.9.0 tools_4.5.0 cachem_1.1.0