R/saveHDF5MultiAssayExperiment.R
HDF5MultiAssayExperiment.Rd
This function takes a MultiAssayExperiment
object and uses the assays
functionality to obtain data matrices out of the experiments. These are
then saved into the .h5
file format. This function relies heavily on
the HDF5Array
package whose installation is required before use.
saveHDF5MultiAssayExpeirment
preserves the classes contained in the
ExperimentList with the exception of matrix
which is
converted to HDF5Matrix
. Internal SummarizedExperiment
assays are
converted to HDF5-backed assays as in
HDF5Array::saveHDF5SummarizedExperiment
. SummarizedExperiment
objects with multiple i
-th assays will have the first assay take
precedence and others assays will be dropped with a warning.
If the first assay in a SummarizedExperiment
contains an array,
the array is preserved in the process of saving and loading the
HDF5-backed MultiAssayExperiment
.
saveHDF5MultiAssayExperiment(
x,
dir = "h5_mae",
prefix = NULL,
replace = FALSE,
chunkdim = NULL,
level = NULL,
as.sparse = NA,
verbose = NA
)
loadHDF5MultiAssayExperiment(dir = "h5_mae", prefix = NULL)
A MultiAssayExperiment object or derivative
The path (as a single string) to the directory where to save the HDF5-based MultiAssayExperiment object or to load it from.
When saving, the directory will be created if it doesn't already exist.
If the directory already exists and no prefix is specified and
replace
is set to TRUE
, then it's replaced with an
empty directory.
An optional prefix to add to the names of the files created
inside dir
. This allows saving more than one object in the same
directory. When the prefix is NULL
, the name of the x
input
MultiAssayExperiment
is used. To avoid the default setting
use an empty character string i.e., ""
. An underscore (_
) is
appended to the prefix when provided; therefore, typical inputs should be
words, e.g., "test".
When no prefix is specified, should a pre-existing directory be replaced with a new empty one? The content of the pre-existing directory will be lost!
The dimensions of the chunks and the compression level to use for writing the assay data to disk.
Passed to the internal calls to writeHDF5Array
.
See ?writeHDF5Array
for more information.
Whether the assay data should be flagged as sparse or not. If set to
NA
(the default), then the specific as.sparse
value to
use for each assay is determined by calling is_sparse()
on them.
Passed to the internal calls to writeHDF5Array
.
See ?writeHDF5Array
for more information and an
IMPORTANT NOTE.
Set to TRUE
to make the function display progress.
In the case of saveHDF5MultiAssayExperiment()
, verbose
is set to NA
by default, in which case verbosity is controlled
by DelayedArray.verbose.block.processing
option. Setting
verbose
to TRUE
or FALSE
overrides the option.
data("miniACC")
testDir <- file.path(tempdir(), "test_mae")
saveHDF5MultiAssayExperiment(
miniACC, dir = testDir, verbose = TRUE, replace = TRUE
)
#> Start writing assay 1/5 to HDF5 file:
#> /tmp/RtmpswreTL/test_mae/miniACC_experiments.h5
#> / Reading and realizing block 1/1 ...
#> OK
#> \ Writing it ...
#> OK
#> Finished writing assay 1/5 to HDF5 file:
#> /tmp/RtmpswreTL/test_mae/miniACC_experiments.h5
#> Start writing assay 2/5 to HDF5 file:
#> /tmp/RtmpswreTL/test_mae/miniACC_experiments.h5
#> / Reading and realizing block 1/1 ...
#> OK
#> \ Writing it ...
#> OK
#> Finished writing assay 2/5 to HDF5 file:
#> /tmp/RtmpswreTL/test_mae/miniACC_experiments.h5
#> Start writing assay 3/5 to HDF5 file:
#> /tmp/RtmpswreTL/test_mae/miniACC_experiments.h5
#> / Reading and realizing block 1/1 ...
#> OK
#> \ Writing it ...
#> OK
#> Finished writing assay 3/5 to HDF5 file:
#> /tmp/RtmpswreTL/test_mae/miniACC_experiments.h5
#> Start writing assay 4/5 to HDF5 file:
#> /tmp/RtmpswreTL/test_mae/miniACC_experiments.h5
#> / Reading and realizing block 1/1 ...
#> OK
#> \ Writing it ...
#> OK
#> Finished writing assay 4/5 to HDF5 file:
#> /tmp/RtmpswreTL/test_mae/miniACC_experiments.h5
#> Start writing assay 5/5 to HDF5 file:
#> /tmp/RtmpswreTL/test_mae/miniACC_experiments.h5
#> / Reading and realizing block 1/1 ...
#> OK
#> \ Writing it ...
#> OK
#> Finished writing assay 5/5 to HDF5 file:
#> /tmp/RtmpswreTL/test_mae/miniACC_experiments.h5
#> Serialize MultiAssayExperiment object to RDS file:
#> /tmp/RtmpswreTL/test_mae/miniACC_mae.rds
## inspect the files in the dir
list.files(testDir)
#> [1] "miniACC_experiments.h5" "miniACC_mae.rds"
loadHDF5MultiAssayExperiment(
dir = testDir
)
#> A MultiAssayExperiment object of 5 listed
#> experiments with user-defined names and respective classes.
#> Containing an ExperimentList class object of length 5:
#> [1] RNASeq2GeneNorm: SummarizedExperiment with 198 rows and 79 columns
#> [2] gistict: SummarizedExperiment with 198 rows and 90 columns
#> [3] RPPAArray: SummarizedExperiment with 33 rows and 46 columns
#> [4] Mutations: HDF5Matrix with 97 rows and 90 columns
#> [5] miRNASeqGene: SummarizedExperiment with 471 rows and 80 columns
#> Functionality:
#> experiments() - obtain the ExperimentList instance
#> colData() - the primary/phenotype DataFrame
#> sampleMap() - the sample coordination DataFrame
#> `$`, `[`, `[[` - extract colData columns, subset, or experiment
#> *Format() - convert into a long or wide DataFrame
#> assays() - convert ExperimentList to a SimpleList of matrices
#> exportClass() - save data to flat files
## remove example files
unlink(testDir, recursive = TRUE)