@lassedochreden

@nilseling
@NilsEling

1 Data and code availability

To follow this tutorial, please visit https://github.com/BodenmillerGroup/demos/tree/main/docs. The compiled .html file of this workshop is hosted at: https://bodenmillergroup.github.io/demos.

To reproduce the analysis, clone the repository:

git clone https://github.com/BodenmillerGroup/demos.git

and open the EuroBioc2023_workshop.Rmd file in the docs folder.

We will need to install the following packages for the workshop:

if (!require("BiocManager", quietly = TRUE))
    install.packages("BiocManager")

BiocManager::install("cytomapper", "cytoviewer")
library(cytomapper)
library(cytoviewer)

2 Introduction

2.1 Highly multiplexed imaging

Highly multiplexed imaging allows simultaneous spatially and single-cell resolved detection of dozens of biological molecules (e.g. proteins) in their native tissue context. As a result, these technologies allow an in-depth analysis of complex systems and diseases such as the tumor microenvironment (Jackson et al. 2020) and type 1 diabetes progression (Damond et al. 2019).

Imaging-based spatial proteomics methods (Moffitt, Lundberg, and Heyn 2022) can be broadly divided into fluorescent cyclic approaches such as tissue-based cyclic immunofluorescence (t-CyCIF) (Lin et al. 2018) and one-step mass-tag based approaches such as multiplexed ion beam imaging (MIBI) (Angelo et al. 2014) and IMC (Giesen et al. 2014).

Across technologies, the acquired data are commonly stored as multi-channel images, where each pixel encodes the abundance of all acquired markers at a specific position in the tissue. Of note, the instructions below will focus on the visualization and exploration of IMC data as an example. However, data from other technologies such as t-CyCIF or MIBI, which produce pixel-level intensities and (optionally) segmentation masks, can be interactively visualized with cytomapper and cytoviewer.

2.1.1 Imaging mass cytometry

IMC, an advancement of CyTOF, combines antibodies tagged with isotopically pure rare earth metals with laser ablation and mass-spectrometry-based detection to produce high-dimensional images (Giesen et al. 2014). It captures the spatial expression of over 40 proteins in parallel at a sub-cellular resolution of 1 μm. Thus, IMC is able to detect cytoplasmic and nuclear localization of proteins.

2.2 Highly multiplexed image analysis

When performing multiplexed image analysis, the user is often faced with a diverse set of computational tools and complex analysis scripts. The main analysis steps, irrespective of the biological question, include 1) Visual inspection of images for quality control, 2) Image pre-processing and segmentation and 3) Single-cell and spatial analysis.

We developed an interoperabale, modularized computational end-to-end workflow (Windhager, Bodenmiller, and Eling 2021) to process and analyze multiplexed imaging data (Figure 1).

The steinbock framework facilitates multi-channel image processing including raw data pre-processing, image segmentation and feature extraction. Data generated by steinbock can be directly read by the imcRtools R/Bioconductor package for data visualization and spatial analysis (Figure 1).

The cytomapper package (Eling et al. 2020) supports image handling and composite as well as segmentation mask visualization, while the newly developed cytoviewer package (Meyer, Eling, and Bodenmiller 2023) allows interactive and easy-to-use image visualization. Both packages integrate into the Bioconductor framework (Gentleman et al. 2004) for single-cell and image analysis leveraging the image handling and analysis strategies from the EBImage Bioconductor package (Pau et al. 2010) and building on commonly used Bioconductor classes such as SingleCellExperiment (Amezquita et al. 2019), SpatialExperiment (Righelli et al. 2022) and CytoImageList (Eling et al. 2020).

Figure 1: Data output after image processing and its representation and visualization in R. During image processing with the steinbock framework, all provided multi-channel images are first segmented to detect individual cells stored as segmentation masks. Based on the segmentation masks, single-cell features including the mean pixel intensity per cell and channel, morphological features and the x-y positions of the cells (regionprops), and spatial cell graphs are extracted. To work with the generated data in R, all multi-channel images are read into a single CytoImageList object, all segmentation masks are read into a single CytoImageList object and all single-cell features of all images are read into a single SpatialExperiment (or SingleCellExperiment) object. Patient and image level features are stored in the elementMetadata slot of each CytoImageList object and the colData slot of the SpatialExperiment object. For image visualization in R, cytomapper (Eling et al. 2020) and cytoviewer (Meyer, Eling, and Bodenmiller 2023) are used.
Figure 1: Data output after image processing and its representation and visualization in R.
During image processing with the steinbock framework, all provided multi-channel images are first segmented to detect individual cells stored as segmentation masks. Based on the segmentation masks, single-cell features including the mean pixel intensity per cell and channel, morphological features and the x-y positions of the cells (regionprops), and spatial cell graphs are extracted. To work with the generated data in R, all multi-channel images are read into a single CytoImageList object, all segmentation masks are read into a single CytoImageList object and all single-cell features of all images are read into a single SpatialExperiment (or SingleCellExperiment) object. Patient and image level features are stored in the elementMetadata slot of each CytoImageList object and the colData slot of the SpatialExperiment object. For image visualization in R, cytomapper (Eling et al. 2020) and cytoviewer (Meyer, Eling, and Bodenmiller 2023) are used.

2.3 Workshop outline

In this workshop, we showcase two related R/Bioconductor packages, cytomapper (Eling et al. 2020) and cytoviewer (Meyer, Eling, and Bodenmiller 2023), that contain user-friendly functions to visualize and explore the multiplexed read-outs and cell-level information obtained by highly multiplexed imaging data.

Outline:

  1. Download and exploration of example IMC data
  2. cytomapper
  3. cytoviewer
  4. Future developments

The visualization approaches presented here were taken from the IMC data analysis book. The book provides much more detailed information on most relevant IMC data analysis steps.

3 Download and exploration of example data

3.1 Download example dataset

Here, we will use an example IMC cancer dataset to showcase the functionality of cytomapper and cytoviewer. This dataset was generated as part of the Integrated iMMUnoprofiling of large adaptive CANcer patient cohort project (immucan.eu) and includes IMC data for 4 cancer patients diagnosed with different tumor types (head and neck cancer, breast cancer, lung cancer and colorectal cancer).

The data input objects were processed as outlined in the IMC data analysis book and can be downloaded from https://zenodo.org/record/7647079.

# Download
download.file("https://zenodo.org/record/7647079/files/spe.rds",
              destfile = "spe.rds")
download.file("https://zenodo.org/record/7647079/files/images.rds",
              destfile = "images.rds")
download.file("https://zenodo.org/record/7647079/files/masks.rds",
              destfile = "masks.rds")

3.2 Read data into R

3.2.1 Images

The image data was stored as a CytoImageList object containing the spillover corrected multi-channel images (n=14). Each image contains 40 channels and each channel represents the pixel-intensities of one marker (proteins for IMC). The proteins for this dataset are immuno-oncology related targets including Ecad, CD8a and CD68, which mark tumor, CD8+ T cells and myeloid cells, respectively.

# Load images
images <- readRDS("images.rds")
images
## CytoImageList containing 14 image(s)
## names(14): Patient1_001 Patient1_002 Patient1_003 Patient2_001 Patient2_002 Patient2_003 Patient2_004 Patient3_001 Patient3_002 Patient3_003 Patient4_005 Patient4_006 Patient4_007 Patient4_008 
## Each image contains 40 channel(s)
## channelNames(40): MPO HistoneH3 SMA CD16 CD38 HLADR CD27 CD15 CD45RA CD163 B2M CD20 CD68 Ido1 CD3 LAG3 / LAG33 CD11c PD1 PDGFRb CD7 GrzB PDL1 TCF7 CD45RO FOXP3 ICOS CD8a CarbonicAnhydrase CD33 Ki67 VISTA CD40 CD4 CD14 Ecad CD303 CD206 cleavedPARP DNA1 DNA2

The elementMetadata slot of the CytoImageList object stores sample and patient ID as well as indication information.

mcols(images)
## DataFrame with 14 rows and 3 columns
##                 sample_id  patient_id  indication
##               <character> <character> <character>
## Patient1_001 Patient1_001    Patient1       SCCHN
## Patient1_002 Patient1_002    Patient1       SCCHN
## Patient1_003 Patient1_003    Patient1       SCCHN
## Patient2_001 Patient2_001    Patient2         BCC
## Patient2_002 Patient2_002    Patient2         BCC
## ...                   ...         ...         ...
## Patient3_003 Patient3_003    Patient3       NSCLC
## Patient4_005 Patient4_005    Patient4         CRC
## Patient4_006 Patient4_006    Patient4         CRC
## Patient4_007 Patient4_007    Patient4         CRC
## Patient4_008 Patient4_008    Patient4         CRC

The channelNames slot stores the names of the targeted proteins.

channelNames(images)
##  [1] "MPO"               "HistoneH3"         "SMA"              
##  [4] "CD16"              "CD38"              "HLADR"            
##  [7] "CD27"              "CD15"              "CD45RA"           
## [10] "CD163"             "B2M"               "CD20"             
## [13] "CD68"              "Ido1"              "CD3"              
## [16] "LAG3 / LAG33"      "CD11c"             "PD1"              
## [19] "PDGFRb"            "CD7"               "GrzB"             
## [22] "PDL1"              "TCF7"              "CD45RO"           
## [25] "FOXP3"             "ICOS"              "CD8a"             
## [28] "CarbonicAnhydrase" "CD33"              "Ki67"             
## [31] "VISTA"             "CD40"              "CD4"              
## [34] "CD14"              "Ecad"              "CD303"            
## [37] "CD206"             "cleavedPARP"       "DNA1"             
## [40] "DNA2"

The CytoImageList object inherits from a SimpleList object and therefore supports common subsetting operations:

images$Patient1_001
## Image 
##   colorMode    : Grayscale 
##   storage.mode : double 
##   dim          : 600 600 40 
##   frames.total : 40 
##   frames.render: 40 
## 
## imageData(object)[1:5,1:6,1]
##              [,1]         [,2]         [,3]         [,4] [,5]     [,6]
## [1,] 0.000000e+00 0.000000e+00 3.876023e-16 0.000000e+00    0 1.000000
## [2,] 0.000000e+00 0.000000e+00 3.068843e+00 0.000000e+00    1 1.202448
## [3,] 0.000000e+00 5.877243e-15 0.000000e+00 0.000000e+00    0 0.000000
## [4,] 1.018369e-14 1.683428e+00 0.000000e+00 1.683428e+00    0 0.000000
## [5,] 1.763591e+00 1.000000e+00 0.000000e+00 9.436896e-16    0 0.000000

Here, individual images are stored in memory as EBImage::Image objects; however the cytomapper package also supports storing images on disk.

3.2.2 Segmentation masks

The segmentation masks were again stored as a CytoImageList object containing one mask (n=14) for each image. Segmentation masks are defined as one-channel images containing integer values for cells and zero for background.

# Load masks
masks <- readRDS("masks.rds")
masks
## CytoImageList containing 14 image(s)
## names(14): Patient1_001 Patient1_002 Patient1_003 Patient2_001 Patient2_002 Patient2_003 Patient2_004 Patient3_001 Patient3_002 Patient3_003 Patient4_005 Patient4_006 Patient4_007 Patient4_008 
## Each image contains 1 channel

It is crucial that the image names are stored in the elementMetadata slot for matching between multi-channel images and segmentation masks.

mcols(images)
## DataFrame with 14 rows and 3 columns
##                 sample_id  patient_id  indication
##               <character> <character> <character>
## Patient1_001 Patient1_001    Patient1       SCCHN
## Patient1_002 Patient1_002    Patient1       SCCHN
## Patient1_003 Patient1_003    Patient1       SCCHN
## Patient2_001 Patient2_001    Patient2         BCC
## Patient2_002 Patient2_002    Patient2         BCC
## ...                   ...         ...         ...
## Patient3_003 Patient3_003    Patient3       NSCLC
## Patient4_005 Patient4_005    Patient4         CRC
## Patient4_006 Patient4_006    Patient4         CRC
## Patient4_007 Patient4_007    Patient4         CRC
## Patient4_008 Patient4_008    Patient4         CRC
mcols(masks)
## DataFrame with 14 rows and 3 columns
##                 sample_id  patient_id  indication
##               <character> <character> <character>
## Patient1_001 Patient1_001    Patient1       SCCHN
## Patient1_002 Patient1_002    Patient1       SCCHN
## Patient1_003 Patient1_003    Patient1       SCCHN
## Patient2_001 Patient2_001    Patient2         BCC
## Patient2_002 Patient2_002    Patient2         BCC
## ...                   ...         ...         ...
## Patient3_003 Patient3_003    Patient3       NSCLC
## Patient4_005 Patient4_005    Patient4         CRC
## Patient4_006 Patient4_006    Patient4         CRC
## Patient4_007 Patient4_007    Patient4         CRC
## Patient4_008 Patient4_008    Patient4         CRC

3.2.3 SpatialExperiment object

The single-cell information was stored in a SpatialExperiment object. Each column entry represents a single-cell and each row entry represents a single marker. We quantify protein abundance as the mean pixel intensity per marker and cell.

# Load spe
spe <- readRDS("spe.rds")
spe
## class: SpatialExperiment 
## dim: 40 47794 
## metadata(5): color_vectors cluster_codes SOM_codes delta_area
##   filterSpatialContext
## assays(2): counts exprs
## rownames(40): MPO HistoneH3 ... DNA1 DNA2
## rowData names(16): channel name ... marker_class used_for_clustering
## colnames(47794): Patient1_001_1 Patient1_001_2 ... Patient4_008_2844
##   Patient4_008_2845
## colData names(33): sample_id ObjectNumber ... patch_id distToCells
## reducedDimNames(8): UMAP TSNE ... seurat UMAP_seurat
## mainExpName: NULL
## altExpNames(0):
## spatialCoords names(2) : Pos_X Pos_Y
## imgData names(1): sample_id

It also contained various metadata information in the colData slot generated during the analysis pipeline including patient-level information (such as indication) and cell-level information (such as cell type and cell area).

# Explore colData
unique(spe$patient_id)
## [1] "Patient1" "Patient2" "Patient3" "Patient4"
unique(spe$sample_id)
##  [1] "Patient1_001" "Patient1_002" "Patient1_003" "Patient2_001" "Patient2_002"
##  [6] "Patient2_003" "Patient2_004" "Patient3_001" "Patient3_002" "Patient3_003"
## [11] "Patient4_005" "Patient4_006" "Patient4_007" "Patient4_008"
unique(spe$indication)
## [1] "SCCHN" "BCC"   "NSCLC" "CRC"
unique(spe$celltype)
##  [1] "Tumor"       "Stroma"      "Myeloid"     "Plasma_cell" "Treg"       
##  [6] "CD8"         "CD4"         "undefined"   "BnTcell"     "Bcell"      
## [11] "Neutrophil"
range(spe$area)
## [1]   5 466

4 cytomapper

The cytomapper R/Bioconductor package was developed to support the handling and visualization of multiple multi-channel images and segmentation masks (Eling et al. 2020). The main data object for image handling is the CytoImageList container to store multi-channel images and segmentation masks.

We will randomly select 3 images for visualization purposes.

# Load package 
library(cytomapper)

# Sample images
set.seed(220517)
cur_id <- sample(unique(spe$sample_id), 3)

cur_spe <- spe[,spe$sample_id %in% cur_id]
cur_images <- images[names(images) %in% cur_id]
cur_masks <- masks[names(masks) %in% cur_id]

4.1 Pixel-level visualization

The following section gives examples for visualizing individual channels or multiple channels as pseudo-color composite images. For this the cytomapper package exports the plotPixels function which expects a CytoImageList object storing one or multiple multi-channel images. In the simplest use case, a single channel can be visualized as follows:

plotPixels(cur_images, 
           colour_by = "Ecad",
           bcg = list(Ecad = c(0, 5, 1)))

The plot above shows the tissue expression of the epithelial tumor marker E-cadherin on the 3 selected images. The bcg parameter (default c(0, 1, 1)) stands for “background”, “contrast”, “gamma” and controls these attributes of the image. This parameter takes a named list where each entry specifies these attributes per channel. The first value of the numeric vector will be added to the pixel intensities (background); pixel intensities will be multiplied by the second entry of the vector (contrast); pixel intensities will be exponentiated by the third entry of the vector (gamma). In most cases, it is sufficient to adjust the second (contrast) entry of the vector.

The following example highlights the visualization of 6 markers (maximum allowed number of markers) at once per image. The markers indicate the spatial distribution of tumor cells (E-cadherin), T cells (CD3), B cells (CD20), CD8+ T cells (CD8a), plasma cells (CD38) and proliferating cells (Ki67).

plotPixels(cur_images, 
           colour_by = c("Ecad", "CD3", "CD20", "CD8a", "CD38", "Ki67"),
           bcg = list(Ecad = c(0, 5, 1),
                      CD3 = c(0, 5, 1),
                      CD20 = c(0, 5, 1),
                      CD8a = c(0, 5, 1),
                      CD38 = c(0, 8, 1),
                      Ki67 = c(0, 5, 1)))

4.1.1 Adjusting colors

The default colors for visualization are chosen by the additive RGB (red, green, blue) color model. For six markers the default colors are: red, green, blue, cyan (green + blue), magenta (red + blue), yellow (green + red). These colors are the easiest to distinguish by eye. However, you can select other colors for each channel by setting the colour parameter:

plotPixels(cur_images, 
           colour_by = c("Ecad", "CD3", "CD20"),
           bcg = list(Ecad = c(0, 5, 1),
                      CD3 = c(0, 5, 1),
                      CD20 = c(0, 5, 1)),
           colour = list(Ecad = c("black", "burlywood1"),
                         CD3 = c("black", "cyan2"),
                         CD20 = c("black", "firebrick1")))

The colour parameter takes a named list in which each entry specifies the colors from which a color gradient is constructed via colorRampPalette. These are usually vectors of length 2 in which the first entry is "black" and the second entry specifies the color of choice.

4.1.2 Outlining cells on images

The following section highlights the combined visualization of pixel- and cell-level information at once.

First, we can outline all cells on composite images to visually check segmnetation quality.

plotPixels(image = cur_images,
           mask = cur_masks,
           img_id = "sample_id",
           colour_by = "DNA1",
           bcg = list(DNA1 = c(0, 5, 1)),
           colour = list(DNA1 = c("black", "blue")),
           thick = TRUE)

In the next step, we can color the cell outlines by a metadata entry. By specifying the outline_by parameter, the outlines of cells can now be colored based on their metadata.

The following example first generates a 3-channel composite images displaying the expression of E-cadherin, CD3 and CD20 before coloring the cells’ outlines by their cell phenotype.

plotPixels(image = cur_images,
           mask = cur_masks,
           object = cur_spe, 
           cell_id = "ObjectNumber", 
           img_id = "sample_id",
           colour_by = c("Ecad", "CD3", "CD20"),
           outline_by = "celltype",
           bcg = list(Ecad = c(0, 5, 1),
                      CD3 = c(0, 5, 1),
                      CD20 = c(0, 5, 1)),
           colour = list(celltype = metadata(cur_spe)$color_vectors$celltype),
           thick = TRUE)

Distinguishing individual cell phenotypes is nearly impossible in the images above.

However, the SpatialExperiment object can be subsetted to only contain cells of a single or few phenotypes. This allows the selective visualization of cell outlines on composite images.

Here, we select all CD8+ T cells from the dataset and outline them on a 2-channel composite image displaying the expression of CD3 and CD8a.

CD8 <- cur_spe[,cur_spe$celltype == "CD8"]

plotPixels(image = cur_images,
           mask = cur_masks,
           object = CD8, 
           cell_id = "ObjectNumber", img_id = "sample_id",
           colour_by = c("CD3", "CD8a"),
           outline_by = "celltype",
           bcg = list(CD3 = c(0, 5, 1),
                      CD8a = c(0, 5, 1)),
           colour = list(celltype = c("CD8" = "white")),
           thick = TRUE)

This type of visualization allows the quality control of two things: 1. segmentation quality of individual cell types can be checked and 2. cell phenotyping accuracy can be visually assessed against expected marker expression.

4.2 Cell-level visualization

In the following section, we will show examples on how to visualize single cells as segmentation masks. This type of visualization allows to observe the spatial distribution of cell phenotypes and the visual assessment of morphological features.

4.2.1 Visualizing metadata

The cytomapper package provides the plotCells function that accepts a CytoImageList object containing segmentation masks. These are defined as single channel images where sets of pixels with the same integer ID identify individual cells. This integer ID can be found as an entry in the colData(spe) slot and as pixel information in the segmentation masks. The entry in colData(spe) needs to be specified via the cell_id argument to the plotCells function. In that way, data contained in the SpatialExperiment object can be mapped to segmentation masks. For the current dataset, the cell IDs are stored in colData(spe)$ObjectNumber.

As cell IDs are only unique within a single image, plotCells also requires the img_id argument. This argument specifies the colData(spe) as well as the mcols(masks) entry that stores the unique image name from which each cell was extracted. In the current dataset the unique image names are stored in colData(spe)$sample_id and mcols(masks)$sample_id.

Providing these two entries that allow mapping between the SpatialExperiment object and segmentation masks, we can now color individual cells based on their cell type:

plotCells(cur_masks,
          object = cur_spe, 
          cell_id = "ObjectNumber", 
          img_id = "sample_id",
          colour_by = "celltype")

If only individual cell types should be visualized, the SpatialExperiment object can be subsetted (e.g., to only contain CD8+ T cells). In the following example CD8+ T cells are colored in red and all other cells that are not contained in the dataset are colored in white (as set by the missing_color argument).

CD8 <- cur_spe[,cur_spe$celltype == "CD8"]

plotCells(cur_masks,
          object = CD8, 
          cell_id = "ObjectNumber", 
          img_id = "sample_id",
          colour_by = "celltype",
          colour = list(celltype = c(CD8 = "red")),
          missing_colour = "white")

In terms of visualizing metadata, any entry in the colData(spe) slot can be visualized. The plotCells function automatically detects if the entry is continuous or discrete. In this fashion, we can now visualize the area of each cell:

plotCells(cur_masks,
          object = cur_spe, 
          cell_id = "ObjectNumber", 
          img_id = "sample_id",
          colour_by = "area")

5 cytoviewer

This section introduces the cytoviewer R/Bioconductor package (Meyer, Eling, and Bodenmiller 2023) for interactive multi-channel image visualization. The cytoviewer package builds on top of the cytomapper Bioconductor package and extends the static visualization strategies provided by cytomapper via an interactive Shiny application.

The graphical user interface of cytoviewer allows intuitive navigation and little coding experience is required to use the package.

5.1 Application overview

The cytoviewer interface is broadly divided into Image-level (Composite and Channels) and Cell-level (Masks) visualization (As seen above). It allows users to overlay individual images with segmentation masks, integrates well with SingleCellExperiment and SpatialExperiment objects for metadata visualization and supports rapid Image downloads (Figure 2B).

Figure 2: cytoviewer interface and functionality. (A) The supported functionality (right) of cytoviewer depends on the data inputs (left). To match information between the objects, cell (cell_id) and image (img_id) identifiers can be provided. SCE/SPE = SingleCellExperiment/SpatialExperiment. (B) The graphical user interface of cytoviewer is divided into a body, header, and sidebar. The body of cytoviewer includes the image viewer, which has three tabs: Composite (Image-level), Channels (Image-level), and Mask (Cell-level). Zooming is supported for Composite and Mask tabs. The package version, R session information, help page, and a drop-down menu for image downloads are located in the header. The sidebar menu has controls for sample selection, image visualization, mask visualization, and general settings. Scale bar: 150 µm (C) cytoviewer supports different viewing modes. Top: The “channels” tab of image-level visualization displays individual channels. Shown are Ecad (magenta), CD8a (cyan), and CD68 (yellow) marking tumor cells, CD8+ T cells, and myeloid cells, respectively. Center: The “composite” tab of image-level visualization visualizes image composites combining multiple channels. These composite images can be overlayed with cell outlines, which can be colored by cell-specific metadata. Shown here are cell outlines colored by cell area (continous value) and cell type (categorical value; tumor cells in white). Channel color settings are as follows for all markers: Contrast: 2,5; Brightness: 1; Gamma: 1.2. Bottom: The “mask” tab can be used to visualize segmentation masks that can be colored by cell-specific metadata. Shown here are segmentation masks colored by cell area (continuous) and cell type (categorical; tumor cells in magenta). Scale bars: 150 µm. (D) “Image appearance” controls can be used to add legends or titles and to change the scale bar length for image-level (top) and cell level (bottom) visualization. The cell-level mask plot depicts tumor (magenta), myeloid (yellow), and CD8+ T cells (cyan). Scale bars: 100 µm. Adapted from (Meyer, Eling, and Bodenmiller 2023)
Figure 2: cytoviewer interface and functionality. (A) The supported functionality (right) of cytoviewer depends on the data inputs (left). To match information between the objects, cell (cell_id) and image (img_id) identifiers can be provided. SCE/SPE = SingleCellExperiment/SpatialExperiment. (B) The graphical user interface of cytoviewer is divided into a body, header, and sidebar. The body of cytoviewer includes the image viewer, which has three tabs: Composite (Image-level), Channels (Image-level), and Mask (Cell-level). Zooming is supported for Composite and Mask tabs. The package version, R session information, help page, and a drop-down menu for image downloads are located in the header. The sidebar menu has controls for sample selection, image visualization, mask visualization, and general settings. Scale bar: 150 µm (C) cytoviewer supports different viewing modes. Top: The “channels” tab of image-level visualization displays individual channels. Shown are Ecad (magenta), CD8a (cyan), and CD68 (yellow) marking tumor cells, CD8+ T cells, and myeloid cells, respectively. Center: The “composite” tab of image-level visualization visualizes image composites combining multiple channels. These composite images can be overlayed with cell outlines, which can be colored by cell-specific metadata. Shown here are cell outlines colored by cell area (continous value) and cell type (categorical value; tumor cells in white). Channel color settings are as follows for all markers: Contrast: 2,5; Brightness: 1; Gamma: 1.2. Bottom: The “mask” tab can be used to visualize segmentation masks that can be colored by cell-specific metadata. Shown here are segmentation masks colored by cell area (continuous) and cell type (categorical; tumor cells in magenta). Scale bars: 150 µm. (D) “Image appearance” controls can be used to add legends or titles and to change the scale bar length for image-level (top) and cell level (bottom) visualization. The cell-level mask plot depicts tumor (magenta), myeloid (yellow), and CD8+ T cells (cyan). Scale bars: 100 µm. Adapted from (Meyer, Eling, and Bodenmiller 2023)

5.1.1 Data input format

The cytoviewer package combines objects of SingleCellExperiment, SpatialExperiment and cytomapper::CytoImageList classes (from cytomapper) to visualize image- and cell-level information.

The cytoviewer function call takes up to five arguments (Figure 2A):

Firstly, image refers to a CytoImageList object containing one or multiple multi-channel images where each channel represents the pixel-intensities of one marker (proteins in IMC).

Secondly, mask refers to a CytoImageList object containing one or multiple segmentation masks. Segmentation masks are defined as one-channel images containing integer values, which represent the cell ids or background.

Thirdly, the object entry refers to a SingleCellExperiment or SpatialExperiment class object that contains cell-specific metadata in the colData slot.

Lastly, to match information between the CytoImageList objects and the SingleCellExperiment/SpatialExperiment object, two additional spots can be specified:

  • img_id: a single character indicating the colData (of the SingleCellExperiment/SpatialExperiment object) and elementMetadata (of the CytoImageList object) entry that contains the image identifiers. These image ids have to match between the SingleCellExperiment/ SpatialExperiment object and the CytoImageList objects.

  • cell_id: a single character indicating the colData entry that contains the cell identifiers. These should be integer values corresponding to pixel-values in the segmentation masks.

5.1.2 Data input variations

The functionality of cytoviewer depends on which input objects are user-provided. Below we describe the four use cases in respect to input objects and functionality (Figure 2A).

1. Usage of cytoviewer with images, masks and object

The full functionality of cytoviewer can be leveraged when image, mask and object are provided, which is the main intended use case.

This allows image-level visualization (Composite and Channels), cell-level visualization, overlaying images with segmentation masks as well as metadata visualization.

2. Usage of cytoviewer with images only

If only the image object is specified, image-level visualization (Composite and Channels) is possible.

3. Usage of cytoviewer with images and masks

Image-level visualization (Composite and Channels), overlaying of images with segmentation masks and Cell-level visualization is feasible when image and mask objects are provided.

4. Usage of cytoviewer with masks and object

If mask and object are specified, cell-level visualization as well as metadata visualization is possible.

5.2 Function call

Here as an example, we will call cytoviewer with the image, mask and object data from the cytomapper section to leverage all provided functionality.

This setting allows image-level visualization (Composite and Channels), cell-level visualization, overlaying images with segmentation masks as well as metadata visualization.

For further details, please refer to the ?cytoviewer manual or the Help page within the shiny application.

library(cytoviewer)

# Use cytoviewer with images, masks and object
app <- cytoviewer(image = cur_images, 
                  mask = cur_masks, 
                  object = cur_spe, 
                  img_id = "sample_id", 
                  cell_id = "ObjectNumber")

if (interactive()) {
  
  shiny::runApp(app)

  }

5.3 Application overview

5.4 Image-level visualization

Image visualization control is split into basic and advanced controls (Figure 2C).

Basic controls supports the selection of up to six markers/channels for image display. Each marker has color control settings that allow the user to set contrast, brightness, gamma and select a channel color.

In the advanced controls part, the user can choose to overlay the displayed images with provided segmentation masks. Outline color and mask thickness can be adjusted by the user. Moreover, the masks can be outlined by cell-specific metadata provided in colData slot of the object.

Of note, for categorical and continuous metadata entries the user can choose between discrete colors and continuous color palettes (viridis, inferno, plasma), respectively.

5.5 Cell-level visualization

Cell visualization has basic controls (Figure 2C).

Here, the user can choose to display the provided segmentation masks. If an object is provided, the masks can be colored by cell-specific metadata.

Please note again that for categorical and continuous metadata entries the user can choose between discrete colors and continuous color palettes (viridis, inferno, plasma), respectively.

5.6 General controls

General controls is subdivided into an Image appearance and Image filters part (Figure 2D).

In the Image appearance section, the user can adjust the scale bar length and include legend/image titles, while the Image filters section allows to control pixel-wise interpolation (default) and apply a Gaussian filter.

5.7 Image download

The cytoviewer package supports fast and uncomplicated image downloads. Download controls are part of the Header.

The user can specify a file name, select the image of interest (Composite, Channels, Mask) and the file format (pdf, png). Upon clicking the download button, a pop-window should appear where the user can specify the download location.

6 Future developments

  1. Future developments of both packages include an integration with modern imaging file types such as OME-NGFF and spatialData to enable image and single-cell analysis in a programming language agnostic fashion.

  2. We further envision to integrate interactive visualization of spatial cell graphs into cytoviewer (“Points-level”).

7 Further resources

The IMC data analysis book contains a detailed overview on the presented and other approaches for multiplexed image analysis and visualization.

The steinbock framework provides functionality for image processing.

The ImcSegmentationPipeline provides a GUI-based version of the segmentation pipeline based on Ilastik pixel classification and image segmentation via CellProfiler.

The imcRtools package supports reading single-cell data from segmented images, multi-channel spillover correction, and spatial data analysis.

The imcdatasets R/Bioconductor package contains a number of publically available IMC datasets.

8 Acknowledgments

Nils Eling, Nicolas Damond and Tobias Hoch developed the cytomapper package. Lasse Meyer and Nils Eling developed the cytoviewer package.

We thank Prof. Bernd Bodenmiller for his continued support.

Nils Eling is funded by Marie Sklodowska Curie Actions.

Session Info

## R version 4.3.0 (2023-04-21)
## Platform: x86_64-apple-darwin20 (64-bit)
## Running under: macOS Big Sur 11.7.9
## 
## Matrix products: default
## BLAS:   /Library/Frameworks/R.framework/Versions/4.3-x86_64/Resources/lib/libRblas.0.dylib 
## LAPACK: /Library/Frameworks/R.framework/Versions/4.3-x86_64/Resources/lib/libRlapack.dylib;  LAPACK version 3.11.0
## 
## locale:
## [1] en_US.UTF-8/en_US.UTF-8/en_US.UTF-8/C/en_US.UTF-8/en_US.UTF-8
## 
## time zone: Europe/Zurich
## tzcode source: internal
## 
## attached base packages:
## [1] stats4    stats     graphics  grDevices utils     datasets  methods  
## [8] base     
## 
## other attached packages:
##  [1] cytoviewer_1.0.0            cytomapper_1.12.0          
##  [3] SingleCellExperiment_1.22.0 SummarizedExperiment_1.30.2
##  [5] Biobase_2.60.0              GenomicRanges_1.52.0       
##  [7] GenomeInfoDb_1.36.1         IRanges_2.34.1             
##  [9] S4Vectors_0.38.1            BiocGenerics_0.46.0        
## [11] MatrixGenerics_1.12.3       matrixStats_1.0.0          
## [13] EBImage_4.42.0              BiocStyle_2.28.0           
## 
## loaded via a namespace (and not attached):
##   [1] splines_4.3.0               later_1.3.1                
##   [3] bitops_1.0-7                tibble_3.2.1               
##   [5] R.oo_1.25.0                 svgPanZoom_0.3.4           
##   [7] polyclip_1.10-4             XML_3.99-0.14              
##   [9] lifecycle_1.0.3             rstatix_0.7.2              
##  [11] edgeR_3.42.4                doParallel_1.0.17          
##  [13] lattice_0.21-8              MASS_7.3-60                
##  [15] backports_1.4.1             magrittr_2.0.3             
##  [17] limma_3.56.2                sass_0.4.7                 
##  [19] rmarkdown_2.23              plotrix_3.8-2              
##  [21] jquerylib_0.1.4             yaml_2.3.7                 
##  [23] httpuv_1.6.11               sp_2.0-0                   
##  [25] cowplot_1.1.1               RColorBrewer_1.1-3         
##  [27] ConsensusClusterPlus_1.64.0 multcomp_1.4-25            
##  [29] abind_1.4-5                 zlibbioc_1.46.0            
##  [31] Rtsne_0.16                  purrr_1.0.1                
##  [33] R.utils_2.12.2              RCurl_1.98-1.12            
##  [35] TH.data_1.1-2               sandwich_3.0-2             
##  [37] tweenr_2.0.2                circlize_0.4.15            
##  [39] GenomeInfoDbData_1.2.10     ggrepel_0.9.3              
##  [41] irlba_2.3.5.1               CATALYST_1.24.0            
##  [43] terra_1.7-39                dqrng_0.3.0                
##  [45] svglite_2.1.1               DelayedMatrixStats_1.22.1  
##  [47] codetools_0.2-19            DropletUtils_1.20.0        
##  [49] DelayedArray_0.26.7         scuttle_1.10.2             
##  [51] ggforce_0.4.1               tidyselect_1.2.0           
##  [53] shape_1.4.6                 raster_3.6-23              
##  [55] farver_2.1.1                ScaledMatrix_1.8.1         
##  [57] viridis_0.6.4               jsonlite_1.8.7             
##  [59] BiocNeighbors_1.18.0        GetoptLong_1.0.5           
##  [61] ellipsis_0.3.2              scater_1.28.0              
##  [63] ggridges_0.5.4              survival_3.5-5             
##  [65] iterators_1.0.14            systemfonts_1.0.4          
##  [67] foreach_1.5.2               tools_4.3.0                
##  [69] ggnewscale_0.4.9            Rcpp_1.0.11                
##  [71] glue_1.6.2                  gridExtra_2.3              
##  [73] xfun_0.40                   dplyr_1.1.2                
##  [75] HDF5Array_1.28.1            shinydashboard_0.7.2       
##  [77] withr_2.5.0                 BiocManager_1.30.22        
##  [79] fastmap_1.1.1               rhdf5filters_1.12.1        
##  [81] fansi_1.0.4                 rsvd_1.0.5                 
##  [83] digest_0.6.33               R6_2.5.1                   
##  [85] mime_0.12                   colorspace_2.1-0           
##  [87] gtools_3.9.4                jpeg_0.1-10                
##  [89] R.methodsS3_1.8.2           utf8_1.2.3                 
##  [91] tidyr_1.3.0                 generics_0.1.3             
##  [93] data.table_1.14.8           htmlwidgets_1.6.2          
##  [95] S4Arrays_1.0.5              pkgconfig_2.0.3            
##  [97] gtable_0.3.3                ComplexHeatmap_2.16.0      
##  [99] RProtoBufLib_2.12.1         XVector_0.40.0             
## [101] htmltools_0.5.5             carData_3.0-5              
## [103] bookdown_0.35               fftwtools_0.9-11           
## [105] clue_0.3-64                 scales_1.2.1               
## [107] png_0.1-8                   SpatialExperiment_1.10.0   
## [109] colorRamps_2.3.1            knitr_1.43                 
## [111] rstudioapi_0.15.0           reshape2_1.4.4             
## [113] rjson_0.2.21                zoo_1.8-12                 
## [115] cachem_1.0.8                rhdf5_2.44.0               
## [117] GlobalOptions_0.1.2         stringr_1.5.0              
## [119] shinycssloaders_1.0.0       miniUI_0.1.1.1             
## [121] parallel_4.3.0              vipor_0.4.5                
## [123] pillar_1.9.0                grid_4.3.0                 
## [125] vctrs_0.6.3                 promises_1.2.0.1           
## [127] ggpubr_0.6.0                BiocSingular_1.16.0        
## [129] car_3.1-2                   cytolib_2.12.1             
## [131] beachmat_2.16.0             xtable_1.8-4               
## [133] cluster_2.1.4               archive_1.1.5              
## [135] beeswarm_0.4.0              evaluate_0.21              
## [137] magick_2.7.5                mvtnorm_1.2-2              
## [139] cli_3.6.1                   locfit_1.5-9.8             
## [141] compiler_4.3.0              rlang_1.1.1                
## [143] crayon_1.5.2                ggsignif_0.6.4             
## [145] FlowSOM_2.8.0               plyr_1.8.8                 
## [147] flowCore_2.12.2             ggbeeswarm_0.7.2           
## [149] stringi_1.7.12              viridisLite_0.4.2          
## [151] BiocParallel_1.34.2         nnls_1.4                   
## [153] munsell_0.5.0               tiff_0.1-11                
## [155] colourpicker_1.2.0          Matrix_1.6-0               
## [157] sparseMatrixStats_1.12.2    ggplot2_3.4.2              
## [159] Rhdf5lib_1.22.0             shiny_1.7.4.1              
## [161] highr_0.10                  drc_3.0-1                  
## [163] fontawesome_0.5.1           memoise_2.0.1              
## [165] igraph_1.5.0.1              broom_1.0.5                
## [167] bslib_0.5.0

References

Amezquita, Robert A., Aaron T. L. Lun, Etienne Becht, Vince J. Carey, Lindsay N. Carpp, Ludwig Geistlinger, Federico Marini, et al. 2019. “Orchestrating Single-Cell Analysis with Bioconductor.” Nature Methods 17 (2): 137–45. https://doi.org/10.1038/s41592-019-0654-x.
Angelo, Michael, Sean C. Bendall, Rachel Finck, Matthew B. Hale, Chuck Hitzman, Alexander D. Borowsky, Richard M. Levenson, et al. 2014. “Multiplexed Ion Beam Imaging of Human Breast Tumors.” Nature Medicine 20 (4): 436–42.
Damond, Nicolas, Stefanie Engler, Vito R. T. Zanotelli, Denis Schapiro, Clive H. Wasserfall, Irina Kusmartseva, Harry S. Nick, et al. 2019. “A Map of Human Type 1 Diabetes Progression by Imaging Mass Cytometry.” Cell Metabolism 29 (3): 755–768.e5.
Eling, Nils, Nicolas Damond, Tobias Hoch, and Bernd Bodenmiller. 2020. “Cytomapper: An r/Bioconductor Package for Visualization of Highly Multiplexed Imaging Data.” Bioinformatics 36 (24): 5706--5708. https://doi.org/10.1093/bioinformatics/btaa1061.
Gentleman, Robert C, Vincent J Carey, Douglas M Bates, Ben Bolstad, Marcel Dettling, Sandrine Dudoit, Byron Ellis, et al. 2004. Genome Biology 5 (10): R80. https://doi.org/10.1186/gb-2004-5-10-r80.
Giesen, Charlotte, Hao A. O. Wang, Denis Schapiro, Nevena Zivanovic, Andrea Jacobs, Bodo Hattendorf, Peter J. Schüffler, et al. 2014. “Highly Multiplexed Imaging of Tumor Tissues with Subcellular Resolution by Mass Cytometry.” Nature Methods 11 (4): 417–22.
Jackson, Hartland W., Jana R. Fischer, Vito R. T. Zanotelli, H. Raza Ali, Robert Mechera, Savas D. Soysal, Holger Moch, et al. 2020. “The Single-Cell Pathology Landscape of Breast Cancer.” Nature 578 (7796): 615–20.
Lin, Jia-Ren, Benjamin Izar, Shu Wang, Clarence Yapp, Shaolin Mei, Parin M. Shah, Sandro Santagata, and Peter K. Sorger. 2018. “Highly Multiplexed Immunofluorescence Imaging of Human Tissues and Tumors Using t-CyCIF and Conventional Optical Microscopes.” eLife 7: 1–46.
Meyer, Lasse, Nils Eling, and Bernd Bodenmiller. 2023. “Cytoviewer: An r/Bioconductor Package for Interactive Visualization and Exploration of Highly Multiplexed Imaging Data.”
Moffitt, Jeffrey R., Emma Lundberg, and Holger Heyn. 2022. “The Emerging Landscape of Spatial Profiling Technologies.” Nature Reviews Genetics 23: 741–59.
Pau, Gregoire, Florian Fuchs, Oleg Sklyar, Michael Boutros, and Wolfgang Huber. 2010. “EBImage: An r Package for Image Processing with Applications to Cellular Phenotypes.” Bioinformatics 26 (7): 979–81. https://doi.org/10.1093/bioinformatics/btq046.
Righelli, Dario, Lukas M Weber, Helena L Crowell, Brenda Pardo, Leonardo Collado-Torres, Shila Ghazanfar, Aaron T L Lun, Stephanie C Hicks, and Davide Risso. 2022. SpatialExperiment: Infrastructure for Spatially-Resolved Transcriptomics Data in r Using Bioconductor.” Edited by Valentina Boeva. Bioinformatics 38 (11): 3128–31. https://doi.org/10.1093/bioinformatics/btac299.
Windhager, Jonas, Bernd Bodenmiller, and Nils Eling. 2021. “An End-to-End Workflow for Multiplexed Image Processing and Analysis.” https://doi.org/10.1101/2021.11.12.468357.