This class facilitates the handling of multiple one- or multi-channel images. It inherits from SimpleList setting elementType="Image". Therefore, each slot contains an either one- or multi-dimensional array in form of an Image object.

CytoImageList(
  ...,
  on_disk = FALSE,
  h5FilesPath = NULL,
  BPPARAM = SerialParam()
)

Arguments

...

A list of images (or coercible to a list) or individual images

on_disk

Logical indicating if images in form of HDF5Array objects (as .h5 files) should be stored on disk rather than in memory.

h5FilesPath

path to where the .h5 files for on disk representation are stored. This path needs to be defined when on_disk = TRUE. When files should only temporarily be stored on disk, please set h5FilesPath = getHDF5DumpDir()

BPPARAM

parameters for parallelised processing. This is only recommended for very large images. See MulticoreParam for information on how to use multiple cores for parallelised processing.

Value

A CytoImageList object

Details

Similar to the Image class, the first two dimensions of each entry indicate the spatial dimension of the image. These can be different for each entry. The third dimension indicates the number of channels per Image. Each entry in the CytoImageList class object must contain the same number of channels. Here, each channel represents pixel values indicating measurement intensities or in case of segmentation masks the cells' ID. The CytoImageList class therefore only supports a Grayscale colormode (see colormode) representation of each individual image.

The class further contains an elementMetadata slot that stores image-level meta information. This slot should be accessed using the mcols accessor function.

Restrictions on entry names

The CytoImageList class only supports unique entry names to avoid duplicated images. Names of a CytoImageList object can be get and set via names(x), where x is a CytoImageList object. Furthermore, only named or unnamed CytoImageList objects are allowed. Partially named objects causing empty or NA names return an error.

Coercion

Coercion to and from list, SimpleList and List:

as.list(x), as(x, "SimpleList"), as(x, "SimpleList"):

Coercion from a CytoImageList object x

as(x, "CytoImageList"):

Coercion from a list, SimpleList or List object x to anCytoImageList object

Looping

While lapply and mapply return regular list objects, endoapply and mendoapply return CytoImageList objects.

On disk representation

When setting on_disk = TRUE and specifying the h5FilesPath, images are stored on disk. To convert back to an in-memory CytoImageList object, one can call CytoImageList(on_disk_IL, on_disk = FLASE).

See also

Image, for further image analysis tools.

SimpleList, for basics functions to handle SimpleList objects

?loadImages, for reading images into a CytoImageList object

?"CytoImageList-naming", for setting and getting image and channel names

?"CytoImageList-subsetting", for subsetting and accessor functions

Author

Nils Eling (nils.eling@dqbm.uzh.ch)

Examples

# Creation of CytoImageList
u <- matrix(rbinom(100, 10, 0.5), ncol=10, nrow=10)
v <- matrix(rbinom(100, 10, 0.5), ncol=10, nrow=10)
IL1 <- CytoImageList(image1 = Image(u), image2 = Image(v))

# Coercion
as.list(IL1)
#> $image1
#> Image 
#>   colorMode    : Grayscale 
#>   storage.mode : integer 
#>   dim          : 10 10 
#>   frames.total : 1 
#>   frames.render: 1 
#> 
#> imageData(object)[1:5,1:6]
#>      [,1] [,2] [,3] [,4] [,5] [,6]
#> [1,]    3    7    4    2    6    5
#> [2,]    7    4    6    4    5    4
#> [3,]    5    2    6    4    6    4
#> [4,]    3    4    4    6    6    5
#> [5,]    1    5    8    5    3    5
#> 
#> $image2
#> Image 
#>   colorMode    : Grayscale 
#>   storage.mode : integer 
#>   dim          : 10 10 
#>   frames.total : 1 
#>   frames.render: 1 
#> 
#> imageData(object)[1:5,1:6]
#>      [,1] [,2] [,3] [,4] [,5] [,6]
#> [1,]    6    6    7    3    3    6
#> [2,]    6    7    5    6    9    4
#> [3,]    5    2    4    7    5    6
#> [4,]    4    7    7    7    4    4
#> [5,]    5    9    7    5    7    5
#> 
as(IL1, "SimpleList")
#> List of length 2
#> names(2): image1 image2
as(list(image1 = Image(u), image2 = Image(v)), "CytoImageList")
#> CytoImageList containing 2 image(s)
#> names(2): image1 image2 
#> Each image contains 1 channel

# On disk representation
IL1 <- CytoImageList(image1 = Image(u), image2 = Image(v),
                     on_disk = TRUE, 
                     h5FilesPath = HDF5Array::getHDF5DumpDir())