HDF5MultiAssayExperiment: Save a MultiAssayExperiment class object to HDF5 and Rds...

HDF5MultiAssayExperimentR Documentation

Save a MultiAssayExperiment class object to HDF5 and Rds files

Description

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.

Usage

saveHDF5MultiAssayExperiment(
  x,
  dir = "h5_mae",
  prefix = NULL,
  replace = FALSE,
  chunkdim = NULL,
  level = NULL,
  as.sparse = NA,
  verbose = NA
)

loadHDF5MultiAssayExperiment(dir = "h5_mae", prefix = NULL)

Arguments

x

A MultiAssayExperiment object or derivative

dir

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.

prefix

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".

replace

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!

chunkdim, level

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.

as.sparse

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.

verbose

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.

Value

  • saveHDF5MultiAssayExperiment: saves an Rds and h5 file to a directory from the input MultiAssayExperiment

  • loadHDF5MultiAssayExperiment: a MultiAssayExperiment object loaded from a folder as saved by saveHDF5MultiAssayExperiment

Examples


data("miniACC")

testDir <- file.path(tempdir(), "test_mae")

saveHDF5MultiAssayExperiment(
    miniACC, dir = testDir, verbose = TRUE, replace = TRUE
)

## inspect the files in the dir
list.files(testDir)

loadHDF5MultiAssayExperiment(
    dir = testDir
)

## remove example files
unlink(testDir, recursive = TRUE)


vjcitn/biocMultiAssay documentation built on Nov. 24, 2024, 3:04 a.m.