Skip to contents

Generate local indicators of spatial association

Usage

lisa(
  cells,
  r = NULL,
  imageID = "imageID",
  cellType = "cellType",
  spatialCoords = c("x", "y"),
  cores = 1,
  window = "convex",
  window.length = NULL,
  whichParallel = "imageID",
  sigma = NULL,
  lisaFunc = "K",
  minLambda = 0.05,
  BPPARAM = BiocParallel::SerialParam(),
  Rs = r
)

Arguments

cells

A SingleCellExperiment, SpatialExperiment or data frame that contains at least the variables x and y, giving the coordinates of each cell, imageID and cellType.

r

A vector of the radii that the measures of association should be calculated.

imageID

The column which contains image identifiers.

cellType

The column which contains the cell types.

spatialCoords

The columns which contain the x and y spatial coordinates.

cores

Number of cores to use for parallel processing, or a BiocParallel MulticoreParam or SerialParam object.

window

Should the window around the regions be 'square', 'convex' or 'concave'.

window.length

A tuning parameter for controlling the level of concavity when estimating concave windows.

whichParallel

Should the function use parallization on the imageID or the cellType.

sigma

A numeric variable used for scaling when filting inhomogeneous L-curves.

lisaFunc

Either "K" or "L" curve.

minLambda

Minimum value for density for scaling when fitting inhomogeneous L-curves.

BPPARAM

{DEPRECATED} A BiocParallel MulticoreParam or SerialParam object.

Rs

{DEPRECATED} A vector of the radii that the measures of association should be calculated.

Value

A matrix of LISA curves

Examples

library(spicyR)
library(SingleCellExperiment)
#> Loading required package: SummarizedExperiment
#> Loading required package: MatrixGenerics
#> Loading required package: matrixStats
#> 
#> Attaching package: ‘MatrixGenerics’
#> The following objects are masked from ‘package:matrixStats’:
#> 
#>     colAlls, colAnyNAs, colAnys, colAvgsPerRowSet, colCollapse,
#>     colCounts, colCummaxs, colCummins, colCumprods, colCumsums,
#>     colDiffs, colIQRDiffs, colIQRs, colLogSumExps, colMadDiffs,
#>     colMads, colMaxs, colMeans2, colMedians, colMins, colOrderStats,
#>     colProds, colQuantiles, colRanges, colRanks, colSdDiffs, colSds,
#>     colSums2, colTabulates, colVarDiffs, colVars, colWeightedMads,
#>     colWeightedMeans, colWeightedMedians, colWeightedSds,
#>     colWeightedVars, rowAlls, rowAnyNAs, rowAnys, rowAvgsPerColSet,
#>     rowCollapse, rowCounts, rowCummaxs, rowCummins, rowCumprods,
#>     rowCumsums, rowDiffs, rowIQRDiffs, rowIQRs, rowLogSumExps,
#>     rowMadDiffs, rowMads, rowMaxs, rowMeans2, rowMedians, rowMins,
#>     rowOrderStats, rowProds, rowQuantiles, rowRanges, rowRanks,
#>     rowSdDiffs, rowSds, rowSums2, rowTabulates, rowVarDiffs, rowVars,
#>     rowWeightedMads, rowWeightedMeans, rowWeightedMedians,
#>     rowWeightedSds, rowWeightedVars
#> Loading required package: GenomicRanges
#> Loading required package: stats4
#> Loading required package: BiocGenerics
#> 
#> Attaching package: ‘BiocGenerics’
#> The following objects are masked from ‘package:stats’:
#> 
#>     IQR, mad, sd, var, xtabs
#> The following objects are masked from ‘package:base’:
#> 
#>     Filter, Find, Map, Position, Reduce, anyDuplicated, aperm, append,
#>     as.data.frame, basename, cbind, colnames, dirname, do.call,
#>     duplicated, eval, evalq, get, grep, grepl, intersect, is.unsorted,
#>     lapply, mapply, match, mget, order, paste, pmax, pmax.int, pmin,
#>     pmin.int, rank, rbind, rownames, sapply, saveRDS, setdiff, table,
#>     tapply, union, unique, unsplit, which.max, which.min
#> Loading required package: S4Vectors
#> 
#> Attaching package: ‘S4Vectors’
#> The following object is masked from ‘package:utils’:
#> 
#>     findMatches
#> The following objects are masked from ‘package:base’:
#> 
#>     I, expand.grid, unname
#> Loading required package: IRanges
#> Loading required package: GenomeInfoDb
#> Loading required package: Biobase
#> Welcome to Bioconductor
#> 
#>     Vignettes contain introductory material; view with
#>     'browseVignettes()'. To cite Bioconductor, see
#>     'citation("Biobase")', and for packages 'citation("pkgname")'.
#> 
#> Attaching package: ‘Biobase’
#> The following object is masked from ‘package:MatrixGenerics’:
#> 
#>     rowMedians
#> The following objects are masked from ‘package:matrixStats’:
#> 
#>     anyMissing, rowMedians
# Read in data
isletFile <- system.file("extdata", "isletCells.txt.gz", package = "spicyR")
cells <- read.table(isletFile, header = TRUE)
cellExp <- SingleCellExperiment(
  assay = list(intensities = t(cells[, grepl(names(cells), pattern = "Intensity_")])),
  colData = cells[, !grepl(names(cells), pattern = "Intensity_")]
)

# Cluster cell types
markers <- t(assay(cellExp, "intensities"))
kM <- kmeans(markers, 8)
colData(cellExp)$cluster <- paste("cluster", kM$cluster, sep = "")

# Generate LISA
lisaCurves <- lisa(
  cellExp,
  spatialCoords = c("Location_Center_X", "Location_Center_Y"),
  cellType = "cluster", imageID = "ImageNumber"
)
#> Generating local L-curves.

# Cluster the LISA curves
kM <- kmeans(lisaCurves, 2)