read_steinbock: Reads in single-cell data generated by the steinbock pipeline

View source: R/read_steinbock.R

read_steinbockR Documentation

Reads in single-cell data generated by the steinbock pipeline

Description

Reader function to generate a SpatialExperiment or SingleCellExperiment object from single-cell data obtained by the steinbock pipeline.

Usage

read_steinbock(
  path,
  intensities_folder = "intensities",
  regionprops_folder = "regionprops",
  graphs_folder = "neighbors",
  pattern = NULL,
  extract_cellid_from = "Object",
  extract_coords_from = c("centroid-1", "centroid-0"),
  image_file = "images.csv",
  extract_imagemetadata_from = c("width_px", "height_px"),
  panel_file = "panel.csv",
  extract_names_from = "name",
  return_as = c("spe", "sce"),
  BPPARAM = SerialParam()
)

Arguments

path

full path to the steinbock output folder

intensities_folder

name of the folder containing the intensity measurements per image

regionprops_folder

name of the folder containing the cell-specific morphology and spatial measurements per image. Can be set to NULL to exclude reading in morphology measures.

graphs_folder

name of the folder containing the spatial connectivity graphs per image. Can be set to NULL to exclude reading in graphs.

pattern

regular expression specifying a subset of files that should be read in.

extract_cellid_from

single character indicating which column entry in the intensity files contains the integer cell id.

extract_coords_from

character vector indicating which column entries in the regionprops files contain the x (first entry) and y (second entry) coordinates.

image_file

single character indicating the file name storing meta data per image (can be NULL).

extract_imagemetadata_from

character vector indicating which additional image specific metadata to extract from the image_file. These will be stored in the colData(x) slot as object/cell-specific entries.

panel_file

single character containing the name of the panel file. This can either be inside the steinbock path (recommended) or located somewhere else.

extract_names_from

single character indicating the column of the panel file containing the channel names.

return_as

should the object be returned as SpatialExperiment (return_as = "spe") or SingleCellExperiment (return_as = "sce").

BPPARAM

parameters for parallelised processing.

Value

returns a SpatialExperiment or SingleCellExperiment object markers in rows and cells in columns.

The returned data container

In the case of both containers x, intensity features are stored in the counts(x) slot. Morphological features are stored in the colData(x) slot. The graphs are stored as SelfHits object in the colPair(x, "neighborhood") slot.

In the case of a returned SpatialExperiment object, the cell coordinates are stored in the spatialCoords(x) slot.

In the case of a returned SingleCellExperiment object, the cell coordinates are stored in the colData(x) slot named as Pos_X and Pos_Y.

Author(s)

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

See Also

https://github.com/BodenmillerGroup/steinbock for the pipeline

read_cpout for reading in single-cell data as produced by the ImcSegmentationPipeline

SingleCellExperiment and SpatialExperiment for the constructor functions.

colPair for information on how to work with the cell-cell interaction graphs

bpparam for the parallelised backend

Examples

path <- system.file("extdata/mockData/steinbock", package = "imcRtools")

# Read in as SpatialExperiment object
x <- read_steinbock(path)
x

# Read in as SingleCellExperiment object
x <- read_steinbock(path, return_as = "sce")
x

# Read in a subset of files
x <- read_steinbock(path, pattern = "mockData1")
x

# Only read in intensities
x <- read_steinbock(path, graphs_folder = NULL, regionprops_folder = NULL)
x

# Parallelisation
x <- read_steinbock(path, BPPARAM = BiocParallel::bpparam())


BodenmillerGroup/imcRtools documentation built on Oct. 30, 2024, 12:19 p.m.