dimGeometries: Dimension geometry methods

dimGeometriesR Documentation

Dimension geometry methods

Description

"Dimension geometry" refers to Simple Feature (sf) geometries associated with rows (features, genes) or columns (cells or spots) of the gene count matrix in the SpatialFeatureExperiment object. For each dimension, the number of rows in the sf data frame specifying the geometries must match the size of the dimension of interest. For example, there must be the same number of rows in the sf data frame describing cells as there are cells in the gene count matrix. This page documents getters and setters for the dimension geometries. The getters and setters are implemented in a way similar to those of reducedDims in SingleCellExperiment.

Usage

## S4 method for signature 'SpatialFeatureExperiment'
dimGeometries(x, MARGIN = 2, withDimnames = TRUE)

## S4 replacement method for signature 'SpatialFeatureExperiment'
dimGeometries(x, MARGIN, withDimnames = TRUE, translate = TRUE, ...) <- value

## S4 method for signature 'SpatialFeatureExperiment'
dimGeometryNames(x, MARGIN)

## S4 replacement method for signature 'SpatialFeatureExperiment,numeric,character'
dimGeometryNames(x, MARGIN) <- value

## S4 method for signature 'SpatialFeatureExperiment,missing'
dimGeometry(x, type, MARGIN, sample_id = NULL, withDimnames = TRUE)

## S4 method for signature 'SpatialFeatureExperiment,numeric'
dimGeometry(x, type, MARGIN, sample_id = NULL, withDimnames = TRUE)

## S4 method for signature 'SpatialFeatureExperiment,character'
dimGeometry(x, type, MARGIN, sample_id = NULL, withDimnames = TRUE)

## S4 replacement method for signature 'SpatialFeatureExperiment,missing'
dimGeometry(
  x,
  type,
  MARGIN,
  sample_id = NULL,
  withDimnames = TRUE,
  translate = TRUE,
  ...
) <- value

## S4 replacement method for signature 'SpatialFeatureExperiment,numeric'
dimGeometry(
  x,
  type,
  MARGIN,
  sample_id = NULL,
  withDimnames = TRUE,
  translate = TRUE,
  ...
) <- value

## S4 replacement method for signature 'SpatialFeatureExperiment,character'
dimGeometry(
  x,
  type,
  MARGIN,
  sample_id = NULL,
  withDimnames = TRUE,
  translate = TRUE,
  ...
) <- value

colGeometry(x, type = 1L, sample_id = NULL, withDimnames = TRUE)

colGeometry(
  x,
  type = 1L,
  sample_id = NULL,
  withDimnames = TRUE,
  translate = TRUE
) <- value

colGeometries(x, withDimnames = TRUE)

colGeometries(x, withDimnames = TRUE, translate = TRUE) <- value

colGeometryNames(x)

colGeometryNames(x) <- value

rowGeometry(x, type = 1L, sample_id = NULL, withDimnames = TRUE)

rowGeometry(
  x,
  type = 1L,
  sample_id = NULL,
  withDimnames = TRUE,
  translate = TRUE
) <- value

rowGeometries(x, withDimnames = TRUE)

rowGeometries(x, withDimnames = TRUE, translate = TRUE) <- value

rowGeometryNames(x)

rowGeometryNames(x) <- value

spotPoly(x, sample_id = NULL, withDimnames = TRUE)

spotPoly(x, sample_id = NULL, withDimnames = TRUE, translate = TRUE) <- value

centroids(x, sample_id = NULL, withDimnames = TRUE)

centroids(x, sample_id = NULL, withDimnames = TRUE, translate = TRUE) <- value

ROIPoly(x, sample_id = NULL, withDimnames = TRUE)

ROIPoly(x, sample_id = NULL, withDimnames = TRUE, translate = TRUE) <- value

cellSeg(x, sample_id = NULL, withDimnames = TRUE)

cellSeg(x, sample_id = NULL, withDimnames = TRUE, translate = TRUE) <- value

nucSeg(x, sample_id = NULL, withDimnames = TRUE)

nucSeg(x, sample_id = NULL, withDimnames = TRUE, translate = TRUE) <- value

txSpots(x, withDimnames = TRUE)

txSpots(x, withDimnames = TRUE, translate = TRUE) <- value

Arguments

x

A SpatialFeatureExperiment object.

MARGIN

As in apply. 1 stands for rows and 2 stands for columns.

withDimnames

Logical. If TRUE, then the dimnames (colnames or rownames) of the gene count matrix should correspond to row names of the sf data frames of interest.

translate

Logical. Only used if removeEmptySpace has been run of the SFE object. If that's the case, this argument indicates whether the new value to be assigned to the geometry is in the coordinates prior to removal of empty space so it should be translated to match the new coordinates after removing empty space. Default to TRUE.

...

spatialCoordsNames, spotDiameter, geometryType passed to df2sf. Defaults are the same as in df2sf. For dimGeometries<- only: geometryType can be a character vector of the geometry type of each data frame in the list of the same length as the list if the data frames specify different types of geometries.

value

Value to set. For dimGeometry, must be a sf data frame with the same number of rows as size in the dimension of interest, or an ordinary data frame that can be converted to such a sf data frame (see df2sf). For dimGeometries, must be a list of such sf or ordinary data frames.

type

An integer specifying the index or string specifying the name of the *Geometry to query or replace. If missing, then the first item in the *Geometries will be returned or replaced.

sample_id

Sample ID to get or set geometries.

Details

These are convenience wrappers for getters and setters of special geometries:

colGeometry/ies

dimGeometry/ies with MARGIN = 2, for geometries associated with columns of the gene count matrix (cells/Visium spots/samples).

rowGeometry/ies

dimGeometry/ies with MARGIN = 1, for geometries associated with rows of the gene count matrix (genes/features).

spotPoly

Polygons of spots from technologies such as Visium, ST, and slide-seq, which do not correspond to cells. Centroids of the polygons are stored in spatialCoords of the underlying SpatialExperiment object.

ROIPoly

Polygons of regions of interest (ROIs) from technologies such as laser capture microdissection (LCM) and GeoMX DSP. These should correspond to columns of the gene count matrix.

cellSeg

Cell segmentation polygons. If the columns of the gene count matrix are single cells, then this is stored in colGeometries. Otherwise, this is stored in annotGeometries.

nucSeg

Similar to cellSeg, but for nuclei rather than whole cell.

txSpots

POINT or MULTIPOINT geometries of transcript spots of single molecular resolution technologies, stored in rowGeometries.

Value

Getters for multiple geometries return a named list. Getters for names return a character vector of the names. Getters for single geometries return an sf data frame. Setters return an SFE object.

Examples

library(SFEData)
sfe <- McKellarMuscleData(dataset = "small")

# Get all column geometries as a named list
# Use MARGIN = 1 or rowGeometry/ies for rowGeometries
cgs <- dimGeometries(sfe, MARGIN = 2)
# Or equivalently
cgs <- colGeometries(sfe)

# Set all column geometries with a named list
dimGeometries(sfe, MARGIN = 2) <- cgs
# Or equivalently
colGeometries(sfe) <- cgs

# Get names of column geometries
cgns <- dimGeometryNames(sfe, MARGIN = 2)
cgns <- colGeometryNames(sfe)

# Set column geometry names
dimGeometryNames(sfe, MARGIN = 2) <- cgns
colGeometryNames(sfe) <- cgns

# Get a specific column geometry by name
spots <- dimGeometry(sfe, "spotPoly", MARGIN = 2)
spots <- colGeometry(sfe, "spotPoly")
# Or equivalently, the wrapper specifically for Visium spot polygons,
# for the name "spotPoly"
spots <- spotPoly(sfe)
# Other colGeometry wrappers for specific names:
# ROIPoly (for LCM and GeoMX DSP), cellSeg and nucSeg (for MERFISH; would
# query annotGeometries for Visium)
# rowGeometry wrappers for specific names: txSpots (MERFISH transcript spots)
# By index
spots <- colGeometry(sfe, 1L)

# Multiple samples, only get geometries for one sample
sfe2 <- McKellarMuscleData("small2")
sfe_combined <- cbind(sfe, sfe2)
spots1 <- colGeometry(sfe, "spotPoly", sample_id = "Vis5A")
spots2 <- spotPoly(sfe_combined, sample_id = "sample02")
# Get geometries for multiple samples
spots3 <- spotPoly(sfe_combined, sample_id = c("Vis5A", "sample02"))
# All samples
spots3 <- spotPoly(sfe_combined, sample_id = "all")

# Set specific column geometry by name
colGeometry(sfe, "foobar") <- spots
# Or use wrapper
spotPoly(sfe) <- spots
# Specify sample_id
colGeometry(sfe_combined, "foobar", sample_id = "Vis5A") <- spots1
# Only entries for the specified sample are set.
foobar <- colGeometry(sfe_combined, "foobar", sample_id = "sample02")

pachterlab/SpatialFeatureExperiment documentation built on March 11, 2024, 11:11 p.m.