R/read_steinbock.R
read_steinbock.Rd
Reader function to generate a
SpatialExperiment
or
SingleCellExperiment
object from single-cell data
obtained by the
steinbock
pipeline.
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()
)
full path to the steinbock output folder
name of the folder containing the intensity measurements per image
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.
name of the folder containing the spatial connectivity
graphs per image. Can be set to NULL
to exclude reading in graphs.
regular expression specifying a subset of files that should be read in.
single character indicating which column entry in the intensity files contains the integer cell id.
character vector indicating which column entries in the regionprops files contain the x (first entry) and y (second entry) coordinates.
single character indicating the file name storing meta data
per image (can be NULL
).
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.
single character containing the name of the panel file. This can either be inside the steinbock path (recommended) or located somewhere else.
single character indicating the column of the panel file containing the channel names.
should the object be returned as
SpatialExperiment
(return_as = "spe"
) or
SingleCellExperiment
(return_as = "sce"
).
parameters for parallelised processing.
returns a SpatialExperiment
or SingleCellExperiment
object markers in rows and cells in columns.
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
.
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
path <- system.file("extdata/mockData/steinbock", package = "imcRtools")
# Read in as SpatialExperiment object
x <- read_steinbock(path)
x
#> class: SpatialExperiment
#> dim: 5 404
#> metadata(0):
#> assays(1): counts
#> rownames(5): Ag107 Cytokeratin 5 Laminin YBX1 H3K27Ac
#> rowData names(7): channel name ... cellpose Tube.Number
#> colnames(404): 20210305_NE_mockData1_001_1 20210305_NE_mockData1_001_2
#> ... 20210305_NE_mockData5_003_39 20210305_NE_mockData5_003_40
#> colData names(8): sample_id ObjectNumber ... width_px height_px
#> reducedDimNames(0):
#> mainExpName: NULL
#> altExpNames(0):
#> spatialCoords names(2) : Pos_X Pos_Y
#> imgData names(1): sample_id
# Read in as SingleCellExperiment object
x <- read_steinbock(path, return_as = "sce")
x
#> class: SingleCellExperiment
#> dim: 5 404
#> metadata(0):
#> assays(1): counts
#> rownames(5): Ag107 Cytokeratin 5 Laminin YBX1 H3K27Ac
#> rowData names(7): channel name ... cellpose Tube.Number
#> colnames(404): 20210305_NE_mockData1_001_1 20210305_NE_mockData1_001_2
#> ... 20210305_NE_mockData5_003_39 20210305_NE_mockData5_003_40
#> colData names(10): sample_id ObjectNumber ... width_px height_px
#> reducedDimNames(0):
#> mainExpName: NULL
#> altExpNames(0):
# Read in a subset of files
x <- read_steinbock(path, pattern = "mockData1")
x
#> class: SpatialExperiment
#> dim: 5 63
#> metadata(0):
#> assays(1): counts
#> rownames(5): Ag107 Cytokeratin 5 Laminin YBX1 H3K27Ac
#> rowData names(7): channel name ... cellpose Tube.Number
#> colnames(63): 20210305_NE_mockData1_001_1 20210305_NE_mockData1_001_2
#> ... 20210305_NE_mockData1_003_24 20210305_NE_mockData1_003_25
#> colData names(8): sample_id ObjectNumber ... width_px height_px
#> reducedDimNames(0):
#> mainExpName: NULL
#> altExpNames(0):
#> spatialCoords names(2) : Pos_X Pos_Y
#> imgData names(1): sample_id
# Only read in intensities
x <- read_steinbock(path, graphs_folder = NULL, regionprops_folder = NULL)
x
#> class: SpatialExperiment
#> dim: 5 404
#> metadata(0):
#> assays(1): counts
#> rownames(5): Ag107 Cytokeratin 5 Laminin YBX1 H3K27Ac
#> rowData names(7): channel name ... cellpose Tube.Number
#> colnames(404): 20210305_NE_mockData1_001_1 20210305_NE_mockData1_001_2
#> ... 20210305_NE_mockData5_003_39 20210305_NE_mockData5_003_40
#> colData names(4): sample_id ObjectNumber width_px height_px
#> reducedDimNames(0):
#> mainExpName: NULL
#> altExpNames(0):
#> spatialCoords names(0) :
#> imgData names(1): sample_id
# Parallelisation
x <- read_steinbock(path, BPPARAM = BiocParallel::bpparam())