measureObjects: Compute morphological and intensity features from objects on...

View source: R/measureObjects.R

measureObjectsR Documentation

Compute morphological and intensity features from objects on images.

Description

For each object (e.g. cell) identified by segmentation, the measureObjects function computes intensity features (also referred to as basic features; e.g. mean intensity), shape features (e.g. area), moment features (e.g. position) and haralick features. These features are returned in form of a SingleCellExperiment or SpatialExperiment object.

Usage

measureObjects(
  mask,
  image,
  img_id,
  return_as = c("sce", "spe"),
  feature_types = c("basic", "shape", "moment"),
  basic_feature = "mean",
  basic_quantiles = NULL,
  shape_feature = c("area", "radius.mean"),
  moment_feature = c("cx", "cy", "majoraxis", "eccentricity"),
  haralick_feature = NULL,
  haralick_nbins = 32,
  haralick_scales = c(1, 2),
  BPPARAM = SerialParam()
)

Arguments

mask

a CytoImageList object containing single-channel Image or HDF5Array objects. Segmentation masks must contain integer pixel values where groups of pixels correspond to objects.

image

a CytoImageList object containing single or multi-channel Image or HDF5Array objects, where each channel indicates the measured pixel intensities.

img_id

character specifying the mcols(image) and mcols(mask) entry, in which the image IDs are stored.

return_as

single character specifying the class of the returned object. This is either "sce" to return a SingleCellExperiment (default) or "spe" to return a SpatialExperiment object.

feature_types

character vector or string indicating which features to compute. Needs to contain "basic". Optionally, "shape", "moment" and "haralick" are allowed. Default "basic", "shape" and "moment".

basic_feature

string indicating which intensity measurement per object and channel should be used to populate the counts(x) slot; where x is the returned object. Default "mean" but "sd", "mad" and "q*" allowed. Here, * indicates the computed quantile (see basic_quantiles).

basic_quantiles

numeric vector or single number indicating which quantiles to compute. Default none.

shape_feature

string or character vector specifying which shape features to compute. Default "area" and "radius.mean". Allowed entries are: "area", "perimeter", "radius.mean", "radius.sd", "radius.max", "radius.min".

moment_feature

string or character vector indicating which moment features to compute. Default "cx", "cy", "majoraxis", and "eccentricity". Other allowed features are "theta". Here moment features are only computed on the segmentation mask without incorporating pixel intensities. Therefore, "cx" and "cy" are the x and y coordinates of the cell centroids.

haralick_feature

string or character vector indicating which haralick features to compute. Default none. Allowed are the 13 haralick features: "asm", "con", "cor", "var", "idm", "sav", "sva", "sen", "ent", "dva", "den", "f12", "f13"

haralick_nbins

an integer indicating the number of bins used to compute the haralick matrix. Pixel intensities are binned in haralick_nbins discrete gray levels before computing the haralick matrix.

haralick_scales

an integer vector indicating the number of scales (distance at which to consider neighbouring pixels) to use to compute the haralick features.

BPPARAM

parameters for parallelised processing of images. See MulticoreParam for information on how to use multiple cores for parallelised processing.

Value

A SingleCellExperiment or SpatialExperiment object (see details)

The returned objects

By default, a SingleCellExperiment object is returned. When setting return_as = "spe", the returned object is of class SpatialExperiment. The returned object contains a single assay. This assay contains individual objects in columns and channels in rows. Each entry summarises the intensities per object and channel. This summary statistic is typically the mean intensity per object and channel. However, other summary statistics can be computed. When the mean intensity per object and channel is computed (default), the assay is accessible via counts(sce). Otherwise, the assay needs to be accessed via assay(sce, "counts_*"), where * indicates the argument to basic_feature.

The colData(x) entry is populated by the computed shape, moment and haralick features per object. The prefix of the feature names indicate whether these features correspond to shape (s.), moment (m.) or haralick (h.) features. Default features are the following:

  • s.area - object size in pixels

  • s.radius.mean - mean object radius in pixels

  • m.cx - x centroid position of object

  • m.cy - y centroid position of object

  • m.majoraxis - major axis length in pixels of elliptical fit

  • m.eccentricity - elliptical eccentricity. 1 meaning straight line and 0 meaning circle.

Computing quantiles

Sometimes it can be useful to describe the summarised pixel intensity per object and channel not in terms of the mean but some quantile of the pixel distribution. For example, to compute the median pixel intensity per object and channel, set basic_feature = "q05" and basic_quantiles = 0.5.

Author(s)

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

See Also

computeFeatures, for detailed explanation for the computed features. https://earlglynn.github.io/RNotes/package/EBImage/Haralick-Textural-Features.html for more discussion on the haralick features

Examples

# Standard example
data(pancreasImages)
data(pancreasMasks)

sce <- measureObjects(pancreasMasks, pancreasImages, img_id = "ImageNb")
sce

# Compute only intensity feature
sce <- measureObjects(pancreasMasks, pancreasImages, img_id = "ImageNb",
                        feature_types = "basic")
colData(sce)

# Visualize on segmentation masks
plotCells(pancreasMasks, object = sce, img_id = "ImageNb",
            cell_id = "object_id", colour_by = "PIN")


BodenmillerGroup/cytomapper documentation built on July 1, 2024, 3:45 p.m.