ImageArrayImageArray 0.99.9
The ImageArray package provides a unified, memory-efficient framework for working with both pyramidal and non-pyramidal images in R by leveraging the DelayedArray package. With ImageArray you can store large images on disk (as HDF5 or Zarr), treat them as array-like objects, perform delayed/lazy operations (e.g., rotate, flip, crop) across all pyramid levels, and seamlessly integrate with other R package and workflows that use large bioimages.
You can install ImageArray from Bioconductor with:
if (!requireNamespace("BiocManager", quietly = TRUE)) {
install.packages("BiocManager")
}
BiocManager::install("ImageArray")
Image pyramids are multi-resolution representations that start with a full resolution image followed by a series of down-sampled levels (e.g., half resolutions). These are common in microscopy, digital pathology and large-scale imaging because:
In our context, an ImageArray object typically contain multiple levels (each a resolution level) internally, and most operations are applied lazily (or delayed) across all levels.
library(ImageArray)
library(BiocFileCache)
library(RBioFormats)
library(EBImage)
library(magick)
library(ggplot2)
library(shiny)
We primarily use images via the EBImage package.
# make random EBImage image
f <- system.file("images", "sample.png", package = "EBImage")
# read image with EBImage
img <- readImage(f)
Similar to its counterparts (see writeHDF5Array in HDF5Array),
we use the writeImageArray function to write pyramid images to preferred
locations. Here, n.levels control the amount of levels that the pyramid image
should have.
output_h5 <- tempfile(fileext = ".h5")
imgarray <- writeImageArray(img, output = output_h5, n.levels = 2)
imgarray
## ImageArray Object (x,y)
## Scales (2): (768,512) (384,256)
Each level of a pyramid can be rasterized at any time, and thus plotted.
bfa.raster <- as.raster(imgarray, level = 2)
plot(bfa.raster)
However, the true functionality of ImageArray is to pick pyramid levels in a memory-efficient manner. Packages, applications and workflows that use ImageArray can rasterize images by defining thresholds for pixel dimensions across all levels.
Here, by using the max.pixel.size, we request ImageArray to
return a pyramid level whose both width (X) and height (Y) are lower than,
for example, 400. This approach is particularly useful when visualizing:
with minimal memory footprint.
# visualize
bfa.raster <- as.raster(imgarray, max.pixel.size = 400)
dim(bfa.raster)
## [1] 256 384
Operations such as rotate, flip, flop or negate will be conducted on all levels (without loading the image in the memory) which then can be rasterized and plotted.
imgarray_rotated <- rotate(imgarray, angle = 90)
imgarray_rotated
## ImageArray Object (x,y)
## Scales (2): (512,768) (256,384)
# visualize
bfa.raster <- as.raster(imgarray_rotated, max.pixel.size = 400)
plot(bfa.raster)
You can even crop images, either by using crop or [ methods.
imgarray_cropped <- imgarray[300:500, 200:300]
imgarray_cropped
## ImageArray Object (x,y)
## Scales (2): (201,101) (101,51)
# visualize
bfa.raster <- as.raster(imgarray_cropped, max.pixel.size = 120)
plot(bfa.raster)
# cropping with indices
imgarray_cropped <- crop(imgarray, index = list(300:500, 200:300))
imgarray_cropped
## ImageArray Object (x,y)
## Scales (2): (201,101) (101,51)
Cropping can also be performed with named indices, thus allowing users
to define indices for only the target dimensions. For now,
ImageArray only includes x and y dimensions and ‘c’ for
image channels.
# cropping with named indices
imgarray_cropped <- crop(imgarray, index = list(x=300:500, y=200:300))
imgarray_cropped
## ImageArray Object (x,y)
## Scales (2): (201,101) (101,51)
ImageArray is also compatible with the magick package where we retain the dimensionality of the image similar to EBImage.
# make random EBImage image
f <- system.file("images", "sample.png", package = "EBImage")
# read image with magick
img <- image_read(f)
# create image array
output_h5 <- tempfile(fileext = ".h5")
imgarray <- writeImageArray(img, output = output_h5)
imgarray
## ImageArray Object (c,x,y)
## Scales (2): (3,768,512) (3,384,256)
# plot
bfa.raster <- as.raster(imgarray, level = 2)
plot(bfa.raster)
You can create ImageArray objects from OME-TIFF files (or any file compatible
with Bio-formats, see RBioFormats) with already defined layers.
ome.tiff.file <- system.file("extdata", "xy_12bit__plant.ome.tiff",
package = "ImageArray")
read.metadata(ome.tiff.file)
## ImageMetadata list of length 2
##
## series res sizeX sizeY sizeC sizeZ sizeT total
## 1 1 512 512 1 1 1 1
## 1 2 256 256 1 1 1 1
##
## globalMetadata:List of 1058
## $ xy_12bit__plant.oir inTrigger triggerNo #5 : chr "0"
## $ xy_12bit__plant.oir inTrigger triggerNo #4 : chr "1"
## $ xy_12bit__plant.oir laser name #7 : chr "LD640"
## $ xy_12bit__plant.oir inTrigger triggerNo #3 : chr "0"
## $ xy_12bit__plant.oir laser name #6 : chr "LD594"
## [list output truncated]
# define ImageArray object
imgarray <- createImageArray(ome.tiff.file, series = 1, resolution = 1:2)
imgarray
## ImageArray Object (x,y,c)
## Scales (2): (512,512,1) (256,256,1)
You can again use either max.pixel.size or min.pixel.size to control
the size of the image that loaded into the memory.
bfa.raster <- as.raster(imgarray, max.pixel.size = 300)
plot(bfa.raster)
dim(bfa.raster)
## [1] 256 256
bfa.raster <- as.raster(imgarray, min.pixel.size = 300)
plot(bfa.raster)
dim(bfa.raster)
## [1] 512 512
bfa.raster <- as.raster(imgarray, level = 1)
plot(bfa.raster)
dim(bfa.raster)
## [1] 512 512
The delayed pyramid scheme introduced by ImageArray objects can also
be used to generate scalable plots and interactive Shiny application
for visualizing large images without loading them into the memory.
For this example, we will use a large H&E image used by 10x Genomics, generated after a Xenium in Situ platform run (i.e. postXenium H&E).
image_file <- paste(
"https://cf.10xgenomics.com/samples/xenium/1.0.1",
"Xenium_FFPE_Human_Breast_Cancer_Rep1",
"Xenium_FFPE_Human_Breast_Cancer_Rep1_he_image.ome.tif",
sep = "/")
library(BiocFileCache)
bfc <- BiocFileCache()
image_file <- bfcrpath(bfc, image_file)
## adding rname 'https://cf.10xgenomics.com/samples/xenium/1.0.1/Xenium_FFPE_Human_Breast_Cancer_Rep1/Xenium_FFPE_Human_Breast_Cancer_Rep1_he_image.ome.tif'
Let us create an ImageArray object first with the image above. The OME-TIFF
file includes multiple resolutions that we can ask ImageArray
object to include.
imgarray <- createImageArray(image_file,
series = 1,
resolution = 1:6)
imgarray
## ImageArray Object (x,y,c)
## Scales (6): (30786,24241,3) (15393,12120,3) ... (1924,1515,3)
## (962,757,3)
You can visualize subsets of large images really quickly using crop and
as.raster. Again, here we use max.pixel.size to optimize the level
parsed from the pyramid; that is, as defined below, only resolutions whose
width and height are both under size 800 would be returned, whether it is
the entire image, or a subset.
# crop image
imgarray_sub <- crop(imgarray, index = list(16000:19000, 7000:10000, NULL))
# convert to raster
img_raster <- as.raster(imgarray_sub, max.pixel.size = 800)
dim(img_raster)
## [1] 751 751
# plot with ggplot
ggplot(data.frame(x = 0, y = 0), aes(x, y)) +
coord_fixed(expand = FALSE,
xlim = c(0, dim(img_raster)[2]),
ylim = c(0, dim(img_raster)[1])) +
annotation_raster(img_raster,
0, dim(img_raster)[2],
dim(img_raster)[1], 0, interpolate = FALSE)
Now let us create a shiny application where we can interactively subset and visualize the OME-TIFF image quickly.
The Shiny application lets the user to define image subsets (or slices) before visualizing without loading large images in memory.
This is because as.raster determines the layer with the pixel dimensions of
he specified subset not exceeding the defined threshold (max.pixel.size).
if (interactive()) {
# variables
dimimg <- dim(imgarray)
max.pixel.size <- 800
# Define UI
ui <- fluidPage(
sliderInput("x_slider", label = "X", min = 1,
max = dimimg[1], value = c(16000, 19000)),
sliderInput("y_slider", label = "Y", min = 1,
max = dimimg[2], value = c(7000, 10000)),
mainPanel(
plotOutput(outputId = "scatterPlot")
)
)
# Define server logic
server <- function(input, output) {
output$scatterPlot <- renderPlot({
# crop image
indx <- seq(input$x_slider[1], input$x_slider[2])
indy <- seq(input$y_slider[1], input$y_slider[2])
imgarray_sub <- crop(imgarray, index = list(indx, indy, NULL))
# convert to raster
img_raster <- as.raster(imgarray_sub, max.pixel.size = max.pixel.size)
# plot with ggplot
imgggplot <- ggplot(data.frame(x = 0, y = 0), aes(x, y)) +
coord_fixed(expand = FALSE,
xlim = c(0, dim(img_raster)[2]),
ylim = c(0, dim(img_raster)[1])) +
annotation_raster(img_raster,
0, dim(img_raster)[2],
dim(img_raster)[1], 0, interpolate = FALSE)
imgggplot
})
}
# Run the app
shinyApp(ui = ui, server = server)
}
## R version 4.6.0 RC (2026-04-17 r89917)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.4 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.23-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=en_US.UTF-8
## [9] LC_ADDRESS=en_US.UTF-8 LC_TELEPHONE=en_US.UTF-8
## [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=en_US.UTF-8
##
## time zone: America/New_York
## tzcode source: system (glibc)
##
## attached base packages:
## [1] stats graphics grDevices utils datasets methods base
##
## other attached packages:
## [1] shiny_1.13.0 ggplot2_4.0.3 magick_2.9.1
## [4] RBioFormats_1.11.0 BiocFileCache_3.1.0 dbplyr_2.5.2
## [7] ImageArray_0.99.9 EBImage_4.53.0 BiocStyle_2.39.0
##
## loaded via a namespace (and not attached):
## [1] tidyselect_1.2.1 dplyr_1.2.1 farver_2.1.2
## [4] blob_1.3.0 S7_0.2.2 filelock_1.0.3
## [7] R.utils_2.13.0 bitops_1.0-9 fastmap_1.2.0
## [10] RCurl_1.98-1.18 promises_1.5.0 digest_0.6.39
## [13] mime_0.13 lifecycle_1.0.5 paws.storage_0.9.0
## [16] RSQLite_2.4.6 magrittr_2.0.5 compiler_4.6.0
## [19] rlang_1.2.0 sass_0.4.10 tools_4.6.0
## [22] yaml_2.3.12 knitr_1.51 S4Arrays_1.11.1
## [25] htmlwidgets_1.6.4 bit_4.6.0 curl_7.1.0
## [28] DelayedArray_0.37.1 RColorBrewer_1.1-3 abind_1.4-8
## [31] HDF5Array_1.39.1 purrr_1.2.2 withr_3.0.2
## [34] BiocGenerics_0.57.1 R.oo_1.27.1 grid_4.6.0
## [37] stats4_4.6.0 xtable_1.8-8 Rhdf5lib_1.99.6
## [40] scales_1.4.0 tinytex_0.59 dichromat_2.0-0.1
## [43] cli_3.6.6 rmarkdown_2.31 crayon_1.5.3
## [46] generics_0.1.4 otel_0.2.0 DBI_1.3.0
## [49] cachem_1.1.0 rhdf5_2.55.16 BiocManager_1.30.27
## [52] XVector_0.51.0 tiff_0.1-12 matrixStats_1.5.0
## [55] vctrs_0.7.3 Matrix_1.7-5 jsonlite_2.0.0
## [58] bookdown_0.46 IRanges_2.45.0 fftwtools_0.9-11
## [61] S4Vectors_0.49.2 bit64_4.8.0 jpeg_0.1-11
## [64] h5mread_1.3.3 locfit_1.5-9.12 jquerylib_0.1.4
## [67] glue_1.8.1 ZarrArray_0.99.5 rJava_1.0-18
## [70] gtable_0.3.6 later_1.4.8 Rarr_1.99.44
## [73] tibble_3.3.1 pillar_1.11.1 rappdirs_0.3.4
## [76] htmltools_0.5.9 rhdf5filters_1.23.3 R6_2.6.1
## [79] httr2_1.2.2 evaluate_1.0.5 lattice_0.22-9
## [82] R.methodsS3_1.8.2 png_0.1-9 memoise_2.0.1
## [85] httpuv_1.6.17 paws.common_0.8.9 bslib_0.10.0
## [88] Rcpp_1.1.1-1.1 SparseArray_1.11.13 xfun_0.57
## [91] MatrixGenerics_1.23.0 pkgconfig_2.0.3