
Identify Spatial EcoTypes from Single-cell Spatial Data (A Single Sample)
Source:R/SpatialEcoTyper.R
SpatialEcoTyper.RdThis function identifies spatially distinct cellular ecosystems (SE) from a single sample.
Usage
SpatialEcoTyper(
normdata,
metadata,
outprefix = "SE",
radius = 50,
resolution = 0.5,
nfeatures = 300,
min.cts.per.region = 2,
npcs = 20,
min.cells = 5,
min.features = 10,
iterations = 10,
minibatch = 5000,
ncores = 4,
grid.size = round(radius * 1.4),
filter.region.by.celltypes = NULL,
k = 20,
k.sn = 50,
dropcell = TRUE
)Arguments
- normdata
A matrix representing normalized gene expression data, where rows correspond to genes and columns correspond to cells.
- metadata
A data frame containing metadata associated with each cell. Must include spatial coordinates (e.g., X and Y) as well as cell type annotations. The row names of the
metadatamust match the column names of thenormdata.- outprefix
Character string specifying the prefix for output file names.
- radius
Numeric specifying the radius (in the same units as spatial coordinates) for defining spatial neighborhoods around each cell. Default is 50.
- resolution
Numeric specifying the resolution for Louvain clustering (default: 0.5).
- nfeatures
Integer specifying the number of top variable features (default: 300) used for the analysis.
- min.cts.per.region
Integer specifying the minimum number of cell types required for a spatial neighborhood.
- npcs
Integer specifying the number of principal components (PCs) (default: 20).
- min.cells
Minimum number of cells / spatial-meta-cells (default: 5) expressing a feature/gene.
- min.features
Minimum number of features (default: 10) detected in a cell / spatial-meta-cell.
- iterations
Integer specifying the number of iterations (default: 10) for SNF analysis.
- minibatch
Integer specifying the number of columns to process in each minibatch in the SNF analysis. Default is 5000. This option splits the matrix into smaller chunks (minibatch), thus reducing memory usage.
- ncores
Integer specifying the number of CPU cores to use for parallel processing.
- grid.size
Numeric specifying the grid size for spatial discretization of coordinates. By default, this size is determined from the specified radius using `round(radius * 1.4)`. Increasing the grid.size will downsample spatial neighborhoods and reduce memory usage, but may remove cells that fall between grid bins.
- filter.region.by.celltypes
A character vector specifying the cell types to include in the analysis. Only spatial neighborhoods that contain at least one of the specified cell types will be analyzed. If NULL, all spatial neighborhoods will be included.
- k
Integer specifying the number of spatial nearest neighbors (default: 20) used to construct spatial meta-cells.
- k.sn
Integer specifying the number of nearest neighbors (default: 50) for constructing similarity networks.
- dropcell
Logical. If TRUE, cells that cannot be assigned to any spatial ecotype will be removed from the returned metadata. Default is TRUE.
Value
A list containing two elements:
- obj
A Seurat object constructed from the fused similarity network of spatial neighborhoods.
- metadata
Updated metadata with a new column (`SE`) containing spatial ecotype labels.
Examples
suppressPackageStartupMessages(library(dplyr))
suppressPackageStartupMessages(library(ggplot2))
suppressPackageStartupMessages(library(parallel))
suppressPackageStartupMessages(library(Seurat))
suppressPackageStartupMessages(library(data.table))
suppressPackageStartupMessages(library(R.utils))
library(SpatialEcoTyper)
# Load single-cell gene expression matrix. Rows are genes, columns are cells.
scdata <- fread("https://spatialecotyper.stanford.edu/inc/inc.public.vignettes.php?file=Melanoma1_subset_counts.tsv.gz",
sep = "\t",header = TRUE, data.table = FALSE)
rownames(scdata) <- scdata[, 1] # set genes as row names
scdata <- as.matrix(scdata[, -1])
normdata <- NormalizeData(scdata)
head(normdata[, 1:5])
#> HumanMelanomaPatient1__cell_3655 HumanMelanomaPatient1__cell_3657
#> PDK4 0.000000 4.221184
#> TNFRSF17 0.000000 0.000000
#> ICAM3 0.000000 0.000000
#> FAP 3.552438 0.000000
#> GZMB 0.000000 0.000000
#> TSC2 0.000000 0.000000
#> HumanMelanomaPatient1__cell_3658 HumanMelanomaPatient1__cell_3660
#> PDK4 3.208945 0
#> TNFRSF17 0.000000 0
#> ICAM3 0.000000 0
#> FAP 0.000000 0
#> GZMB 0.000000 0
#> TSC2 0.000000 0
#> HumanMelanomaPatient1__cell_3661
#> PDK4 0
#> TNFRSF17 0
#> ICAM3 0
#> FAP 0
#> GZMB 0
#> TSC2 0
# Load single-cell metadata. Three columns are required, including X, Y, and CellType. Others are optional.
scmeta <- read.table("https://spatialecotyper.stanford.edu/inc/inc.public.vignettes.php?file=Melanoma1_subset_scmeta.tsv",
sep = "\t",header = TRUE, row.names = 1)
scmeta <- scmeta[colnames(scdata), ] # match the cell ids in scdata and scmeta
head(scmeta)
#> X Y CellType CellTypeName
#> HumanMelanomaPatient1__cell_3655 1894.706 -6367.766 Fibroblast Fibroblasts
#> HumanMelanomaPatient1__cell_3657 1942.480 -6369.602 Fibroblast Fibroblasts
#> HumanMelanomaPatient1__cell_3658 1963.007 -6374.026 Fibroblast Fibroblasts
#> HumanMelanomaPatient1__cell_3660 1981.600 -6372.266 Fibroblast Fibroblasts
#> HumanMelanomaPatient1__cell_3661 1742.939 -6374.851 Fibroblast Fibroblasts
#> HumanMelanomaPatient1__cell_3663 1921.683 -6383.309 Fibroblast Fibroblasts
#> Region Dist2Interface
#> HumanMelanomaPatient1__cell_3655 Stroma -883.1752
#> HumanMelanomaPatient1__cell_3657 Stroma -894.8463
#> HumanMelanomaPatient1__cell_3658 Stroma -904.1115
#> HumanMelanomaPatient1__cell_3660 Stroma -907.8909
#> HumanMelanomaPatient1__cell_3661 Stroma -874.2712
#> HumanMelanomaPatient1__cell_3663 Stroma -903.6559
se_results <- SpatialEcoTyper(normdata, scmeta,
outprefix = "Melanoma1_subset",
radius = 50, ncores = 2)
#> 2026-06-02 16:46:30.090529 Remove 88 genes expressed in fewer than 5 cells
#> 2026-06-02 16:46:38.553194 Construct spatial meta cells for each cell type
#> Total spatial neighborhoods: 2315; total spatial meta cells: 9853
#> 2026-06-02 16:46:39.242726 PCA for each cell type
#> 2026-06-02 16:46:40.932983 Construct cell-type-specific similarity network
#> Removing 182 spatial neighborhoods with fewer than 2 distinct cell types.
#> 2026-06-02 16:46:41.674361 Similarity network fusion
#> 2026-06-02 16:46:41.675671 Normalize networks ...
#> 2026-06-02 16:46:41.776789 Calculate the local transition matrix ...
#> 2026-06-02 16:46:52.402473 Perform the diffusion ...
#> 2026-06-02 16:46:52.404357 Iteration: 1
#> 2026-06-02 16:46:57.972885 Iteration: 2
#> 2026-06-02 16:47:13.942033 Iteration: 3
#> 2026-06-02 16:47:35.454652 Iteration: 4
#> 2026-06-02 16:47:56.332459 Iteration: 5
#> 2026-06-02 16:48:15.005129 Iteration: 6
#> 2026-06-02 16:48:30.708274 Iteration: 7
#> 2026-06-02 16:48:50.098438 Iteration: 8
#> 2026-06-02 16:49:09.46884 Iteration: 9
#> 2026-06-02 16:49:24.502733 Iteration: 10
#> 2026-06-02 16:49:48.42179 Create Seurat object and perform clustering analysis
#> Warning: Data is of class dgeMatrix. Coercing to dgCMatrix.
#> Warning: The default method for RunUMAP has changed from calling Python UMAP via reticulate to the R-native UWOT using the cosine metric
#> To use Python UMAP via reticulate, set umap.method to 'umap-learn' and metric to 'correlation'
#> This message will be shown once per session
#> Error in AnnotateCells(metadata, obj, dropcell = dropcell): The provided Seurat object does not appear to be a valid SpatialEcoTyper result. It must include 'Spot.X' and 'Spot.Y' in obj@meta.data.
# Extract the Seurat object and updated single-cell metadata
obj <- se_results$obj # A Seurat object
#> Error in eval(expr, envir, enclos): object 'se_results' not found
obj
#> Error in eval(expr, envir, enclos): object 'obj' not found
scmeta <- se_results$metadata %>% arrange(SE) # Updated single-cell meta data, with SE annotation added
#> Error in eval(expr, envir, enclos): object 'se_results' not found
head(scmeta)
#> X Y CellType CellTypeName
#> HumanMelanomaPatient1__cell_3655 1894.706 -6367.766 Fibroblast Fibroblasts
#> HumanMelanomaPatient1__cell_3657 1942.480 -6369.602 Fibroblast Fibroblasts
#> HumanMelanomaPatient1__cell_3658 1963.007 -6374.026 Fibroblast Fibroblasts
#> HumanMelanomaPatient1__cell_3660 1981.600 -6372.266 Fibroblast Fibroblasts
#> HumanMelanomaPatient1__cell_3661 1742.939 -6374.851 Fibroblast Fibroblasts
#> HumanMelanomaPatient1__cell_3663 1921.683 -6383.309 Fibroblast Fibroblasts
#> Region Dist2Interface
#> HumanMelanomaPatient1__cell_3655 Stroma -883.1752
#> HumanMelanomaPatient1__cell_3657 Stroma -894.8463
#> HumanMelanomaPatient1__cell_3658 Stroma -904.1115
#> HumanMelanomaPatient1__cell_3660 Stroma -907.8909
#> HumanMelanomaPatient1__cell_3661 Stroma -874.2712
#> HumanMelanomaPatient1__cell_3663 Stroma -903.6559