annotGeometries: Annotation geometry methods

annotGeometriesR Documentation

Annotation geometry methods

Description

"Annotation geometry" refers to Simple Feature (sf) geometries NOT associated with rows (features, genes) or columns (cells or spots) of the gene count matrix in the SpatialFeatureExperiment object. So there can be any number of rows in the sf data frame specifying the geometry. Examples of such geometries are tissue boundaries, pathologist annotation of histological regions, and objects not characterized by columns of the gene count matrix (e.g. nuclei segmentation in a Visium dataset where the columns are Visium spots). This page documents getters and setters for the annotation geometries. Internally, annotation geometries are stored in int_metadata.

Usage

## S4 method for signature 'SpatialFeatureExperiment'
annotGeometries(x)

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

## S4 method for signature 'SpatialFeatureExperiment'
annotGeometryNames(x)

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

## S4 method for signature 'SpatialFeatureExperiment,missing'
annotGeometry(x, type, sample_id = NULL)

## S4 method for signature 'SpatialFeatureExperiment,numeric'
annotGeometry(x, type, sample_id = NULL)

## S4 method for signature 'SpatialFeatureExperiment,character'
annotGeometry(x, type, sample_id = NULL)

## S4 replacement method for signature 'SpatialFeatureExperiment,missing'
annotGeometry(x, type, sample_id = NULL) <- value

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

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

tissueBoundary(x, sample_id = NULL)

tissueBoundary(x, sample_id = NULL, translate = TRUE, ...) <- value

Arguments

x

A SpatialFeatureExperiment object.

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 annotGeometry, must be a sf data frame, or an ordinary data frame that can be converted to a sf data frame (see df2sf). For annotGeometries, must be a list of such sf or ordinary data frames. There must be a column sample_id to indicate the sample the geometries are for, and the sample_id must also appear in colData.

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

Wrapper for getter and setter of special geometry:

tisseuBoundary

Boundary of the tissue of interest, including holes. This is usually of geometry type MULTIPOLYGON, though geometries in annotGeometries can have any type supported by sf.

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

# Example dataset
library(SFEData)
sfe_small <- McKellarMuscleData(dataset = "small")

# Get all annotation geometries, returning a named list
annotGeometries(sfe_small)

# Set all annotation geometries, in a named list
toy <- readRDS(system.file("extdata/sfe_toy.rds",
    package = "SpatialFeatureExperiment"
))
ag <- readRDS(system.file("extdata/ag.rds",
    package = "SpatialFeatureExperiment"
))
annotGeometries(toy) <- list(hull = ag)

# Get names of annotation geometries
annotGeometryNames(sfe_small)

# Set names of annotation geometries
annotGeometryNames(toy) <- "foo"

# Get a specific annotation geometry by name
# sample_id is optional when there is only one sample present
nuclei <- annotGeometry(sfe_small, type = "nuclei", sample_id = "Vis5A")

# Get a specific annotation geometry by index
tb <- annotGeometry(sfe_small, type = 1L)

# Set a specific annotation geometry
annotGeometry(sfe_small, type = "nuclei2") <- nuclei

# Special convenience function for tissue boundaries
# Getter
tb <- tissueBoundary(sfe_small, sample_id = "Vis5A")
# Setter
tissueBoundary(sfe_small, sample_id = "Vis5A") <- tb

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