| Title: | A package for the analysis of anatomical tissue structures in spatial omics data |
|---|---|
| Description: | sosta (Spatial Omics STructure Analysis) is a package for analyzing spatial omics data to explore tissue organization at the anatomical structure level. It reconstructs anatomically relevant structures based on molecular features or cell types. It further calculates a range of metrics at the structure level to quantitatively describe tissue architecture. The package is designed to integrate with other packages for the analysis of spatial omics data. |
| Authors: | Samuel Gunz [aut, cre] (ORCID: <https://orcid.org/0000-0002-8909-0932>), Mark D. Robinson [aut, fnd] |
| Maintainer: | Samuel Gunz <[email protected]> |
| License: | GPL (>= 3) + file LICENSE |
| Version: | 1.5.1 |
| Built: | 2026-06-02 18:49:50 UTC |
| Source: | https://github.com/bioc/sosta |
data.frame to ppp objectAssumes that the data.frame is the output of .SPE2df(). Column order is important!
.df2ppp(df, xName, yName, marks = NULL).df2ppp(df, xName, yName, marks = NULL)
df |
data.frame; with x, y coordinates, image, and categorical mark information. |
xName |
character; column name of x coordinate |
yName |
character; column name of y coordinate |
marks |
character; column name of the mark variable |
ppp; object of type ppp
data(sostaSPE) df <- .SPE2df(sostaSPE, marks = "cellType", imageCol = "imageName") ppp <- .df2ppp(df, xName = "x", yName = "y", marks = "cellType")data(sostaSPE) df <- .SPE2df(sostaSPE, marks = "cellType", imageCol = "imageName") ppp <- .df2ppp(df, xName = "x", yName = "y", marks = "cellType")
Function to estimate the intensity image of a point pattern
.intensityImage(ppp, markSelect = NULL, bndw = NULL, dim).intensityImage(ppp, markSelect = NULL, bndw = NULL, dim)
ppp |
point pattern object of class |
markSelect |
character; name of mark that is to be selected for the reconstruction |
bndw |
bandwidth of kernel density estimator |
dim |
numeric; x dimension of the final reconstruction. |
list; list with the intensity image and the bandwidth and dimension parameters
Function to estimate the intensity threshold for the reconstruction of spatial structures
.intensityThreshold(densityImage, steps = 250, minRange = 0.15).intensityThreshold(densityImage, steps = 250, minRange = 0.15)
densityImage |
real-valued pixel image; output from the function |
steps |
numeric; value used to filter the density estimates, where only
densities greater than the maximum value divided by |
minRange |
numeric; value used to filter minimal threshold in percent of total range. Should range between (0,1]. |
numeric; estimated threshold
SpatialExperiment object to a data frameFunction to convert SpatialExperiment object to a data frame
.SPE2df(spe, imageCol = NULL, marks = NULL, colNames = FALSE).SPE2df(spe, imageCol = NULL, marks = NULL, colNames = FALSE)
spe |
SpatialExperiment; a object of class |
imageCol |
character; name of a column in |
marks |
character; name of column in |
colNames |
logical; extract |
data.frame with x, y coordinates, image, and categorical mark information
data(sostaSPE) .SPE2df(sostaSPE, marks = "cellType", imageCol = "imageName") |> head()data(sostaSPE) .SPE2df(sostaSPE, marks = "cellType", imageCol = "imageName") |> head()
This function assigns each spatial coordinate in a SpatialExperiment object (spe) to the first intersecting structure from a given set of spatial structures.
assingCellsToStructures( spe, allStructs, imageCol = NULL, uniqueId = "structID", nCores = 1 )assingCellsToStructures( spe, allStructs, imageCol = NULL, uniqueId = "structID", nCores = 1 )
spe |
SpatialExperiment; An object of class |
allStructs |
sf; A simple feature collection (sf object) representing spatial structures. Must contain a column which contains a unique identifier for each structure. Default = |
imageCol |
character; The column name in |
uniqueId |
character; The column name in the simple feature collection for which to compute the assignment. |
nCores |
integer; The number of cores to use for parallel processing (default is 1). |
A named list with structure assignments for each spatial point in spe. Points that do not overlap with any structure are assigned NA. Names correspond to colnames of the SpatialExperiment input object.
library("SpatialExperiment") data("sostaSPE") allStructs <- reconstructShapeDensitySPE(sostaSPE, marks = "cellType", imageCol = "imageName", markSelect = "A", bndw = 3.5, thres = 0.045 ) # The function `assingCellsToStructures` needs colnames so we create them here colnames(sostaSPE) <- paste0("cell_", c(1:dim(sostaSPE)[2])) res <- assingCellsToStructures( spe = sostaSPE, allStructs = allStructs, imageCol = "imageName" ) # Assign the structure assignment in the order of the columns in the `SpatialExperiment` object colData(sostaSPE)$structAssign <- res[colnames(sostaSPE)] if (require("ggplot2")) { cbind( colData(sostaSPE[, sostaSPE[["imageName"]] == "image1"]), spatialCoords(sostaSPE[, sostaSPE[["imageName"]] == "image1"]) ) |> as.data.frame() |> ggplot(aes(x = x, y = y, color = structAssign)) + geom_point(size = 0.25) + coord_equal() }library("SpatialExperiment") data("sostaSPE") allStructs <- reconstructShapeDensitySPE(sostaSPE, marks = "cellType", imageCol = "imageName", markSelect = "A", bndw = 3.5, thres = 0.045 ) # The function `assingCellsToStructures` needs colnames so we create them here colnames(sostaSPE) <- paste0("cell_", c(1:dim(sostaSPE)[2])) res <- assingCellsToStructures( spe = sostaSPE, allStructs = allStructs, imageCol = "imageName" ) # Assign the structure assignment in the order of the columns in the `SpatialExperiment` object colData(sostaSPE)$structAssign <- res[colnames(sostaSPE)] if (require("ggplot2")) { cbind( colData(sostaSPE[, sostaSPE[["imageName"]] == "image1"]), spatialCoords(sostaSPE[, sostaSPE[["imageName"]] == "image1"]) ) |> as.data.frame() |> ggplot(aes(x = x, y = y, color = structAssign)) + geom_point(size = 0.25) + coord_equal() }
Converts a binary matrix to an sf polygon
binaryImageToSF(binaryMatrix, xmin, xmax, ymin, ymax)binaryImageToSF(binaryMatrix, xmin, xmax, ymin, ymax)
binaryMatrix |
matrix; binary matrix |
xmin |
integer; minimum x coordinate of the coordinate system |
xmax |
integer; maximum x coordinate of the coordinate system |
ymin |
integer; minimum y coordinate of the coordinate system |
ymax |
integer; maximum y coordinate of the coordinate system |
sf object
matrixR <- matrix(c( 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 ), nrow = 9, byrow = TRUE) polyR <- binaryImageToSF(matrixR, xmin = 0, xmax = 1, ymin = 0, ymax = 1) plot(polyR)matrixR <- matrix(c( 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 ), nrow = 9, byrow = TRUE) polyR <- binaryImageToSF(matrixR, xmin = 0, xmax = 1, ymin = 0, ymax = 1) plot(polyR)
Calculate the proportion of each cell type within spatial structures
cellTypeProportions(spe, structColumn, cellTypeColumn, nCores = 1)cellTypeProportions(spe, structColumn, cellTypeColumn, nCores = 1)
spe |
SpatialExperiment object |
structColumn |
character; name of the |
cellTypeColumn |
character; name of the |
nCores |
integer; The number of cores to use for parallel processing (default is 1). |
A data frame where rows correspond to unique structures and columns correspond to cell types, containing the proportion of each cell type within each structure.
library("SpatialExperiment") data("sostaSPE") allStructs <- reconstructShapeDensitySPE(sostaSPE, marks = "cellType", imageCol = "imageName", markSelect = "A", bndw = 3.5, thres = 0.045 ) # The function `assingCellsToStructures` needs colnames so we create them here colnames(sostaSPE) <- paste0("cell_", c(1:dim(sostaSPE)[2])) # Assign the structure assignment in the order of the columns in the `SpatialExperiment` object colData(sostaSPE)$structAssign <- assingCellsToStructures( spe = sostaSPE, allStructs = allStructs, imageCol = "imageName" )[colnames(sostaSPE)] cellTypeProportions(sostaSPE, "structAssign", "cellType")library("SpatialExperiment") data("sostaSPE") allStructs <- reconstructShapeDensitySPE(sostaSPE, marks = "cellType", imageCol = "imageName", markSelect = "A", bndw = 3.5, thres = 0.045 ) # The function `assingCellsToStructures` needs colnames so we create them here colnames(sostaSPE) <- paste0("cell_", c(1:dim(sostaSPE)[2])) # Assign the structure assignment in the order of the columns in the `SpatialExperiment` object colData(sostaSPE)$structAssign <- assingCellsToStructures( spe = sostaSPE, allStructs = allStructs, imageCol = "imageName" )[colnames(sostaSPE)] cellTypeProportions(sostaSPE, "structAssign", "cellType")
This function generates a spatial point pattern with different types of points (A, B, C) distributed over the simulated tissue structure.
createPointPatternTissue( tissueImage, intA, intB, noiseA = 0.005, intCInA, intCInB )createPointPatternTissue( tissueImage, intA, intB, noiseA = 0.005, intCInA, intCInB )
tissueImage |
Matrix; A binary matrix representing the simulated tissue. |
intA |
Numeric; Intensity of type "A" points (points per unit area) on tissue regions. |
intB |
Numeric; Intensity of type "B" points (points per unit area) on non-tissue regions. |
noiseA |
Numeric; Intensity of type "A" points (points per unit area) on non-tissue regions. |
intCInA |
Numeric; Intensity of type "C" points placed in extended regions around tissue. |
intCInB |
Numeric; Intensity of type "C" points placed within tissue. |
A ppp object representing the spatial point pattern.
tissueImage <- simulateTissueBlobs(128, 100, 7) createPointPatternTissue(tissueImage, 0.01, 0.01, 0.005, 0.005, 0.005)tissueImage <- simulateTissueBlobs(128, 100, 7) createPointPatternTissue(tissueImage, 0.01, 0.01, 0.005, 0.005, 0.005)
Estimate reconstruction parameters from a set of images
estimateReconstructionParametersSPE( spe, marks, imageCol = NULL, markSelect = NULL, nImages = NULL, fun = "bw.diggle", dim = 500, nCores = 1, plotHist = TRUE )estimateReconstructionParametersSPE( spe, marks, imageCol = NULL, markSelect = NULL, nImages = NULL, fun = "bw.diggle", dim = 500, nCores = 1, plotHist = TRUE )
spe |
SpatialExperiment; a object of class |
marks |
character; name of column in |
imageCol |
character; name of a column in |
markSelect |
character; name of mark that is to be selected for the reconstruction |
nImages |
integer; number of images for the estimation. Will be randomly sampled |
fun |
character; function to estimate the kernel density. Default bw.diggle. |
dim |
numeric; x dimension of the final reconstruction. A lower resolution speed up computation but lead to less exact reconstruction. Default = 500 |
nCores |
numeric; number of cores for parallel processing using |
plotHist |
logical; if histogram of estimated densities and thresholds
should be plotted (only if |
tibble; tibble with estimated intensities and thresholds
data("sostaSPE") estimateReconstructionParametersSPE(sostaSPE, marks = "cellType", imageCol = "imageName", markSelect = "A", plotHist = TRUE )data("sostaSPE") estimateReconstructionParametersSPE(sostaSPE, marks = "cellType", imageCol = "imageName", markSelect = "A", plotHist = TRUE )
Estimate the intensity threshold for the reconstruction of spatial structures
findIntensityThreshold(ppp, markSelect = NULL, bndw = NULL, dim, steps = 250)findIntensityThreshold(ppp, markSelect = NULL, bndw = NULL, dim, steps = 250)
ppp |
point pattern object of class |
markSelect |
character; name of mark that is to be selected for the reconstruction |
bndw |
numeric; bandwith of the sigma parameter in the density estimation,
if no value is given the bandwith is estimated using cross validation with
the |
dim |
numeric; x dimension of the final reconstruction. |
steps |
numeric; value used to filter the density estimates, where only
densities greater than the maximum value divided by |
numeric; estimated intensity threshold
data(sostaSPE) ppp <- SPE2ppp(sostaSPE, marks = "cellType", imageCol = "imageName", imageId = "image1") findIntensityThreshold(ppp, markSelect = "A", dim = 250)data(sostaSPE) ppp <- SPE2ppp(sostaSPE, marks = "cellType", imageCol = "imageName", imageId = "image1") findIntensityThreshold(ppp, markSelect = "A", dim = 250)
Function to get the dimension based on dim of y axis
getDimXY(ppp, ydim)getDimXY(ppp, ydim)
ppp |
point pattern object of class |
ydim |
dimension of y axis |
vector; vector with x and y dimension
data(sostaSPE) pp <- SPE2ppp(sostaSPE, marks = "cellType", imageCol = "imageName", imageId = "image1" ) getDimXY(pp, 500)data(sostaSPE) pp <- SPE2ppp(sostaSPE, marks = "cellType", imageCol = "imageName", imageId = "image1" ) getDimXY(pp, 500)
Calculate mean shape metrics of a set of polygons
meanShapeMetrics(totalShapeMetricMatrix)meanShapeMetrics(totalShapeMetricMatrix)
totalShapeMetricMatrix |
matrix of shape metrics |
matrix; matrix of mean shape metrics
data(sostaSPE) struct <- reconstructShapeDensityImage(sostaSPE, marks = "cellType", imageCol = "imageName", imageId = "image1", markSelect = "A", dim = 500 ) shapeMetrics <- totalShapeMetrics(struct) meanShapeMetrics(shapeMetrics)data(sostaSPE) struct <- reconstructShapeDensityImage(sostaSPE, marks = "cellType", imageCol = "imageName", imageId = "image1", markSelect = "A", dim = 500 ) shapeMetrics <- totalShapeMetrics(struct) meanShapeMetrics(shapeMetrics)
Compute minimum boundary distances for each cell within its corresponding image structures
minBoundaryDistances( spe, imageCol, structColumn = NULL, allStructs, nCores = 1 )minBoundaryDistances( spe, imageCol, structColumn = NULL, allStructs, nCores = 1 )
spe |
SpatialExperiment object |
imageCol |
character; name of the |
structColumn |
character; name of the |
allStructs |
sf object; contains spatial structures with corresponding image names |
nCores |
integer; The number of cores to use for parallel processing (default is 1). |
A named list containing the minimum distances between cells and structure boundaries,
values within structures have negative values. Names correspond to colnames of the SpatialExperiment input object.
library("SpatialExperiment") data("sostaSPE") allStructs <- reconstructShapeDensitySPE(sostaSPE, marks = "cellType", imageCol = "imageName", markSelect = "A", bndw = 3.5, thres = 0.045 ) # The function `assingCellsToStructures` needs colnames so we create them here colnames(sostaSPE) <- paste0("cell_", c(1:dim(sostaSPE)[2])) # Assign the structure assignment in the order of the columns in the `SpatialExperiment` object colData(sostaSPE)$structAssign <- assingCellsToStructures( spe = sostaSPE, allStructs = allStructs, imageCol = "imageName" )[colnames(sostaSPE)] res <- minBoundaryDistances( spe = sostaSPE, imageCol = "imageName", structColumn = "structAssign", allStructs = allStructs ) colData(sostaSPE)$minDist <- res[colnames(sostaSPE)] if (require("ggplot2")) { cbind(colData(sostaSPE), spatialCoords(sostaSPE)) |> as.data.frame() |> ggplot(aes(x = x, y = y, color = minDist)) + geom_point(size = 0.25) + scale_colour_gradient2() + geom_sf(data = allStructs, fill = NA, inherit.aes = FALSE) + facet_wrap(~imageName) }library("SpatialExperiment") data("sostaSPE") allStructs <- reconstructShapeDensitySPE(sostaSPE, marks = "cellType", imageCol = "imageName", markSelect = "A", bndw = 3.5, thres = 0.045 ) # The function `assingCellsToStructures` needs colnames so we create them here colnames(sostaSPE) <- paste0("cell_", c(1:dim(sostaSPE)[2])) # Assign the structure assignment in the order of the columns in the `SpatialExperiment` object colData(sostaSPE)$structAssign <- assingCellsToStructures( spe = sostaSPE, allStructs = allStructs, imageCol = "imageName" )[colnames(sostaSPE)] res <- minBoundaryDistances( spe = sostaSPE, imageCol = "imageName", structColumn = "structAssign", allStructs = allStructs ) colData(sostaSPE)$minDist <- res[colnames(sostaSPE)] if (require("ggplot2")) { cbind(colData(sostaSPE), spatialCoords(sostaSPE)) |> as.data.frame() |> ggplot(aes(x = x, y = y, color = minDist)) + geom_point(size = 0.25) + scale_colour_gradient2() + geom_sf(data = allStructs, fill = NA, inherit.aes = FALSE) + facet_wrap(~imageName) }
Compute minimum distances from each cell types to structure boundaries per structure
minCellTypeStructDist( spe, allStructs, structID = "structID", cellTypeColumn, imageCol, nCores = 1 )minCellTypeStructDist( spe, allStructs, structID = "structID", cellTypeColumn, imageCol, nCores = 1 )
spe |
SpatialExperiment object |
allStructs |
sf object; contains spatial structures with corresponding image names |
structID |
character; name of the column in |
cellTypeColumn |
character; name of the |
imageCol |
character; name of the |
nCores |
integer; Number of cores for parallel processing (default = 1) |
A data frame where rows are structure IDs and cols are cell types with minimum distances
data("sostaSPE") allStructs <- reconstructShapeDensitySPE(sostaSPE, marks = "cellType", imageCol = "imageName", markSelect = "A", bndw = 3.5, thres = 0.045 ) minCellTypeStructDist(sostaSPE, allStructs, cellTypeColumn = "cellType", imageCol = "imageName")data("sostaSPE") allStructs <- reconstructShapeDensitySPE(sostaSPE, marks = "cellType", imageCol = "imageName", markSelect = "A", bndw = 3.5, thres = 0.045 ) minCellTypeStructDist(sostaSPE, allStructs, cellTypeColumn = "cellType", imageCol = "imageName")
Function to normalize coordinates between zero and one while keep scaling
normalizeCoordinates(coords)normalizeCoordinates(coords)
coords |
matrix; matrix with coordinates |
matrix; coordinates scaled between 0 and 1
matrixR <- matrix(c( 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 ), nrow = 9, byrow = TRUE) coords <- xyCoordinates(matrixR) normalizeCoordinates(coords)matrixR <- matrix(c( 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 ), nrow = 9, byrow = TRUE) coords <- xyCoordinates(matrixR) normalizeCoordinates(coords)
This function estimates the density of a spatial point pattern (ppp),
thresholds the density to create a binary image, and then converts it
to a valid sf object (polygons).
reconstructShapeDensity( ppp, markSelect = NULL, bndw = NULL, thres = NULL, complement = FALSE, dim )reconstructShapeDensity( ppp, markSelect = NULL, bndw = NULL, thres = NULL, complement = FALSE, dim )
ppp |
point pattern object of class |
markSelect |
character; name of mark that is to be selected for the reconstruction |
bndw |
bandwidth of kernel density estimator |
thres |
intensity threshold for the reconstruction |
complement |
logical; reconstruct everything but the mark of interest, default = |
dim |
numeric; x dimension of the final reconstruction. |
sf object of class POLYGON
data("sostaSPE") ppp <- SPE2ppp(sostaSPE, marks = "cellType", imageCol = "imageName", imageId = "image1") thres <- findIntensityThreshold(ppp, markSelect = "A", dim = 500) struct <- reconstructShapeDensity(ppp, markSelect = "A", thres = thres, dim = 500) plot(struct)data("sostaSPE") ppp <- SPE2ppp(sostaSPE, marks = "cellType", imageCol = "imageName", imageId = "image1") thres <- findIntensityThreshold(ppp, markSelect = "A", dim = 500) struct <- reconstructShapeDensity(ppp, markSelect = "A", thres = thres, dim = 500) plot(struct)
Reconstruct structure from spe object with given image id
reconstructShapeDensityImage( spe, marks, imageCol = NULL, imageId = NULL, markSelect, dim = 500, bndw = NULL, thres = NULL, complement = FALSE )reconstructShapeDensityImage( spe, marks, imageCol = NULL, imageId = NULL, markSelect, dim = 500, bndw = NULL, thres = NULL, complement = FALSE )
spe |
SpatialExperiment; a object of class |
marks |
character; name of column in |
imageCol |
character; name of a column in |
imageId |
character; image id, must be present in imageCol |
markSelect |
character; name of mark that is to be selected for the reconstruction |
dim |
numeric; x dimension of the final reconstruction. A lower resolution speed up computation but lead to less exact reconstruction. Default = 500 |
bndw |
numeric; smoothing bandwidth in the density estimation,
corresponds to the |
thres |
numeric; intensity threshold for the reconstruction; if NULL the threshold is set as the mean between the mode of the pixel intensity distributions |
complement |
logical; reconstruct everything but the mark of interest, default = |
sf object of class POLYGON
data("sostaSPE") struct <- reconstructShapeDensityImage(sostaSPE, marks = "cellType", imageCol = "imageName", imageId = "image1", markSelect = "A", dim = 500 ) plot(struct)data("sostaSPE") struct <- reconstructShapeDensityImage(sostaSPE, marks = "cellType", imageCol = "imageName", imageId = "image1", markSelect = "A", dim = 500 ) plot(struct)
Reconstruct structure from spatial experiment object per image id
reconstructShapeDensitySPE( spe, marks, imageCol = NULL, markSelect, dim = 500, bndw = NULL, thres = NULL, complement = FALSE, nCores = 1 )reconstructShapeDensitySPE( spe, marks, imageCol = NULL, markSelect, dim = 500, bndw = NULL, thres = NULL, complement = FALSE, nCores = 1 )
spe |
SpatialExperiment; a object of class |
marks |
character; name of column in |
imageCol |
character; name of a column in |
markSelect |
character; name of mark that is to be selected for the reconstruction |
dim |
numeric; x dimension of the final reconstruction. A lower resolution speed up computation but lead to less exact reconstruction. Default = 500 |
bndw |
numeric; bandwidth of the sigma parameter in the density estimation,
if no value is given the bandwidth is estimated using cross validation with
the |
thres |
numeric; intensity threshold for the reconstruction; if NULL the threshold is set as the mean between the mode of the pixel intensity distributions estimated for each image individual |
complement |
logical; reconstruct everything but the mark of interest, default = |
nCores |
numeric; number of cores for parallel processing using
|
simple feature collection
data("sostaSPE") allStructs <- reconstructShapeDensitySPE(sostaSPE, marks = "cellType", imageCol = "imageName", markSelect = "A", bndw = 3.5, thres = 0.005 ) allStructsdata("sostaSPE") allStructs <- reconstructShapeDensitySPE(sostaSPE, marks = "cellType", imageCol = "imageName", markSelect = "A", bndw = 3.5, thres = 0.005 ) allStructs
This function plots the intensity of a point pattern image and displays a histogram of the intensity values. Note that intensities less than largest intensity value divided by 250 are not displayed in the histogram.
shapeIntensityImage( spe, marks, imageCol = NULL, imageId = NULL, markSelect, bndw = NULL, dim = 500 )shapeIntensityImage( spe, marks, imageCol = NULL, imageId = NULL, markSelect, bndw = NULL, dim = 500 )
spe |
SpatialExperiment; a object of class |
marks |
character; name of column in |
imageCol |
character; name of a column in |
imageId |
character; image id, must be present in imageCol |
markSelect |
character; name of mark that is to be selected for the reconstruction |
bndw |
numeric; smoothing bandwidth in the density estimation,
corresponds to the |
dim |
numeric; x dimension of the final reconstruction. A lower resolution speeds up computation but lead to less exact reconstruction. Default = 500 |
ggplot object with intensity image and histogram
data("sostaSPE") shapeIntensityImage(sostaSPE, marks = "cellType", imageCol = "imageName", imageId = "image1", markSelect = "A" )data("sostaSPE") shapeIntensityImage(sostaSPE, marks = "cellType", imageCol = "imageName", imageId = "image1", markSelect = "A" )
Calculate a set of shape metrics of a single polygon
shapeMetrics(sfPoly)shapeMetrics(sfPoly)
sfPoly |
POLYGON of class sfc |
For multiple polyogns or a MULTIPOLYGON object use the function totalShapeMetrics
list; list of shape metrics
matrix_R <- matrix(c( 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 ), nrow = 9, byrow = TRUE) polyR <- binaryImageToSF(matrix_R, xmin = 0, xmax = 1, ymin = 0, ymax = 1) shapeMetrics(polyR)matrix_R <- matrix(c( 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 ), nrow = 9, byrow = TRUE) polyR <- binaryImageToSF(matrix_R, xmin = 0, xmax = 1, ymin = 0, ymax = 1) shapeMetrics(polyR)
This function generates a simulated tissue-like structure using a Gaussian blur technique.
simulateTissueBlobs(size, seedNumber, clumpSize)simulateTissueBlobs(size, seedNumber, clumpSize)
size |
Integer; The size (width and height) of the simulated tissue image. |
seedNumber |
Integer; The number of random seed points used to generate tissue blobs. |
clumpSize |
Numeric; The standard deviation (sigma) of the Gaussian blur applied to generate tissue clumps. |
A binary matrix representing the simulated tissue structure.
tissueImage <- simulateTissueBlobs(128, 100, 7) image(tissueImage)tissueImage <- simulateTissueBlobs(128, 100, 7) image(tissueImage)
This dataset contains a simulated SpatialExperiment object (sostaSPE) representing
three tissue images, each with a corresponding spatial point pattern. The point patterns
contain different cell types (A, B, and C), distributed according to simulated
tissue structures.
sostaSPEsostaSPE
A SpatialExperiment object with the following structure:
Numeric; x-coordinate of each point (cell location).
Numeric; y-coordinate of each point (cell location).
Factor; Cell type assigned to each point (A, B, or C).
Factor; Identifier for the tissue image (image1, image2, or image3).
The dataset was generated as follows:
Three tissue images were simulated using simulateTissueBlobs().
Spatial point patterns were created for each tissue using createPointPatternTissue().
The point pattern data was converted into a SpatialExperiment object with spatial coordinates.
Function to convert spatialCoords to an sf object
spatialCoords2SF(spe)spatialCoords2SF(spe)
spe |
SpatialExperiment; a object of class |
sf; Simple feature collection of geometry type POINT
data(sostaSPE) speSel <- sostaSPE[, sostaSPE[["imageName"]] == "image1"] spatialCoords2SF(speSel)data(sostaSPE) speSel <- sostaSPE[, sostaSPE[["imageName"]] == "image1"] spatialCoords2SF(speSel)
SpatialExperiment object to a ppp objectFunction to convert spatial coordinates of a SpatialExperiment object to a ppp object
SPE2ppp(spe, marks = NULL, imageCol = NULL, imageId = NULL)SPE2ppp(spe, marks = NULL, imageCol = NULL, imageId = NULL)
spe |
SpatialExperiment; a object of class |
marks |
character; name of column in |
imageCol |
character; name of a column in |
imageId |
character; image id, must be present in imageCol |
ppp; object of type ppp
data(sostaSPE) SPE2ppp(sostaSPE, marks = "cellType", imageCol = "imageName", imageId = "image1" )data(sostaSPE) SPE2ppp(sostaSPE, marks = "cellType", imageCol = "imageName", imageId = "image1" )
Calculate curvature of sf object
stCalculateCurvature(sfPoly, smoothness = 5)stCalculateCurvature(sfPoly, smoothness = 5)
sfPoly |
|
smoothness |
list; curvature measures |
list; list of curvatures values
https://stackoverflow.com/questions/62250151/calculate-curvature-of-a-closed-object-in-r
matrixR <- matrix(c( 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 ), nrow = 9, byrow = TRUE) polyR <- binaryImageToSF(matrixR, xmin = 0, xmax = 1, ymin = 0, ymax = 1) stCalculateCurvature(polyR)matrixR <- matrix(c( 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 ), nrow = 9, byrow = TRUE) polyR <- binaryImageToSF(matrixR, xmin = 0, xmax = 1, ymin = 0, ymax = 1) stCalculateCurvature(polyR)
Calculate curl of a polygon
stCalculateShapeCurl(sfPoly)stCalculateShapeCurl(sfPoly)
sfPoly |
|
numeric; the curl of the polygon
matrixR <- matrix(c( 1, 1, 1, 1, 1, 0, 1, 1, 0, 0, 1, 1, 1, 1, 0, 0, 1, 1, 1, 1, 1, 1, 1, 0, 1, 1, 0, 1, 1, 0, 1, 1, 0, 0, 1, 1, 1, 1, 0, 0, 1, 1 ), nrow = 7, byrow = TRUE) polyR <- binaryImageToSF(matrixR, xmin = 0, xmax = 1, ymin = 0, ymax = 1) stCalculateShapeCurl(polyR)matrixR <- matrix(c( 1, 1, 1, 1, 1, 0, 1, 1, 0, 0, 1, 1, 1, 1, 0, 0, 1, 1, 1, 1, 1, 1, 1, 0, 1, 1, 0, 1, 1, 0, 1, 1, 0, 0, 1, 1, 1, 1, 0, 0, 1, 1 ), nrow = 7, byrow = TRUE) polyR <- binaryImageToSF(matrixR, xmin = 0, xmax = 1, ymin = 0, ymax = 1) stCalculateShapeCurl(polyR)
Calculate the length of feature axes of a single sf polygon
stFeatureAxes(sfPoly)stFeatureAxes(sfPoly)
sfPoly |
|
list; list containing the major and minor axis lengths
matrixR <- matrix(c( 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 ), nrow = 9, byrow = TRUE) polyR <- binaryImageToSF(matrixR, xmin = 0, xmax = 1, ymin = 0, ymax = 1) stFeatureAxes(polyR)matrixR <- matrix(c( 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 ), nrow = 9, byrow = TRUE) polyR <- binaryImageToSF(matrixR, xmin = 0, xmax = 1, ymin = 0, ymax = 1) stFeatureAxes(polyR)
Calculate a set of shape metrics of a set of polygons
totalShapeMetrics(sfInput)totalShapeMetrics(sfInput)
sfInput |
|
Calculate a set of shape metrics of a set of polygons.
The function calculates all metrics that are implemented in the function
shapeMetrics
matrix; matrix of shape metrics
data(sostaSPE) struct <- reconstructShapeDensityImage(sostaSPE, marks = "cellType", imageCol = "imageName", imageId = "image1", markSelect = "A", dim = 500 ) totalShapeMetrics(struct)data(sostaSPE) struct <- reconstructShapeDensityImage(sostaSPE, marks = "cellType", imageCol = "imageName", imageId = "image1", markSelect = "A", dim = 500 ) totalShapeMetrics(struct)
Function to extract x y coordinates from binary image
xyCoordinates(inputMatrix)xyCoordinates(inputMatrix)
inputMatrix |
a binary matrix |
matrix; matrix with x,y coordinates of the cell of the input matrix
matrixR <- matrix(c( 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 ), nrow = 9, byrow = TRUE) xyCoordinates(matrixR)matrixR <- matrix(c( 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 ), nrow = 9, byrow = TRUE) xyCoordinates(matrixR)