MultiAssayExperiment
s to
artifacts and back againThe alabaster.mae
package implements methods to save MultiAssayExperiment
objects to file artifacts and load them back into R. Check out the
alabaster.base
for more details on the motivation and concepts of the
alabaster framework.
Let’s create a mildly complicated MAE containing RNA-seq and ChIP-seq data with partial overlaps:
library(SummarizedExperiment)
rna.counts <- matrix(rpois(60, 10), ncol=6)
colnames(rna.counts) <- c("disease1", "disease2", "disease3", "control1", "control2", "control3")
rownames(rna.counts) <- c("ENSMUSG00000000001", "ENSMUSG00000000003", "ENSMUSG00000000028",
"ENSMUSG00000000031", "ENSMUSG00000000037", "ENSMUSG00000000049", "ENSMUSG00000000056",
"ENSMUSG00000000058", "ENSMUSG00000000078", "ENSMUSG00000000085")
rna.se <- SummarizedExperiment(list(counts=rna.counts))
colData(rna.se)$disease <- rep(c("disease", "control"), each=3)
chip.counts <- matrix(rpois(100, 10), ncol=4)
colnames(chip.counts) <- c("disease1", "disease2", "control1", "control3")
chip.peaks <- GRanges("chr1", IRanges(1:25*100+1, 1:25*100+100))
chip.se <- SummarizedExperiment(list(counts=chip.counts), rowRanges=chip.peaks)
library(MultiAssayExperiment)
mapping <- DataFrame(
assay = rep(c("rnaseq", "chipseq"), c(ncol(rna.se), ncol(chip.se))), # experiment name
primary = c(colnames(rna.se), colnames(chip.se)), # sample identifiers
colname = c(colnames(rna.se), colnames(chip.se)) # column names inside each experiment
)
mae <- MultiAssayExperiment(list(rnaseq=rna.se, chipseq=chip.se), sampleMap=mapping)
We can use saveObject()
to save it inside a staging
directory:
library(alabaster.mae)
tmp <- tempfile()
meta <- saveObject(mae, tmp)
list.files(tmp, recursive=TRUE)
## [1] "OBJECT"
## [2] "experiments/0/OBJECT"
## [3] "experiments/0/assays/0/OBJECT"
## [4] "experiments/0/assays/0/array.h5"
## [5] "experiments/0/assays/names.json"
## [6] "experiments/0/column_data/OBJECT"
## [7] "experiments/0/column_data/basic_columns.h5"
## [8] "experiments/0/row_data/OBJECT"
## [9] "experiments/0/row_data/basic_columns.h5"
## [10] "experiments/1/OBJECT"
## [11] "experiments/1/assays/0/OBJECT"
## [12] "experiments/1/assays/0/array.h5"
## [13] "experiments/1/assays/names.json"
## [14] "experiments/1/column_data/OBJECT"
## [15] "experiments/1/column_data/basic_columns.h5"
## [16] "experiments/1/row_ranges/OBJECT"
## [17] "experiments/1/row_ranges/ranges.h5"
## [18] "experiments/1/row_ranges/sequence_information/OBJECT"
## [19] "experiments/1/row_ranges/sequence_information/info.h5"
## [20] "experiments/names.json"
## [21] "sample_data/OBJECT"
## [22] "sample_data/basic_columns.h5"
## [23] "sample_map.h5"
We can then load it back into the session with
readObject()
.
## [1] "MultiAssayExperiment"
## attr(,"package")
## [1] "MultiAssayExperiment"
## R version 4.4.1 (2024-06-14)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.1 LTS
##
## Matrix products: default
## BLAS: /usr/lib/x86_64-linux-gnu/openblas-pthread/libblas.so.3
## LAPACK: /usr/lib/x86_64-linux-gnu/openblas-pthread/libopenblasp-r0.3.26.so; LAPACK version 3.12.0
##
## locale:
## [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
## [3] LC_TIME=en_US.UTF-8 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: Etc/UTC
## tzcode source: system (glibc)
##
## attached base packages:
## [1] stats4 stats graphics grDevices utils datasets methods
## [8] base
##
## other attached packages:
## [1] alabaster.mae_1.7.0 alabaster.base_1.5.10
## [3] MultiAssayExperiment_1.31.5 SummarizedExperiment_1.35.5
## [5] Biobase_2.65.1 GenomicRanges_1.57.2
## [7] GenomeInfoDb_1.41.2 IRanges_2.39.2
## [9] S4Vectors_0.43.2 BiocGenerics_0.51.3
## [11] MatrixGenerics_1.17.1 matrixStats_1.4.1
## [13] BiocStyle_2.33.1
##
## loaded via a namespace (and not attached):
## [1] sass_0.4.9 SparseArray_1.5.45 lattice_0.22-6
## [4] alabaster.se_1.5.4 digest_0.6.37 evaluate_1.0.1
## [7] grid_4.4.1 fastmap_1.2.0 jsonlite_1.8.9
## [10] Matrix_1.7-1 alabaster.schemas_1.5.0 BiocManager_1.30.25
## [13] httr_1.4.7 UCSC.utils_1.1.0 HDF5Array_1.33.8
## [16] jquerylib_0.1.4 abind_1.4-8 cli_3.6.3
## [19] rlang_1.1.4 crayon_1.5.3 XVector_0.45.0
## [22] cachem_1.1.0 DelayedArray_0.31.14 yaml_2.3.10
## [25] S4Arrays_1.5.11 tools_4.4.1 Rhdf5lib_1.27.0
## [28] GenomeInfoDbData_1.2.13 alabaster.ranges_1.5.2 alabaster.matrix_1.5.10
## [31] buildtools_1.0.0 R6_2.5.1 lifecycle_1.0.4
## [34] rhdf5_2.49.0 zlibbioc_1.51.2 bslib_0.8.0
## [37] Rcpp_1.0.13 xfun_0.48 sys_3.4.3
## [40] knitr_1.48 rhdf5filters_1.17.0 htmltools_0.5.8.1
## [43] rmarkdown_2.28 maketools_1.3.1 compiler_4.4.1