SPOTlight: Deconvolution of mixture using single-cell data

View source: R/SPOTlight.R

SPOTlightR Documentation

Deconvolution of mixture using single-cell data

Description

This is the backbone function which takes in single cell expression data to deconvolute spatial transcriptomics spots.

Usage

SPOTlight(
  x,
  y,
  groups = NULL,
  mgs,
  n_top = NULL,
  gene_id = "gene",
  group_id = "cluster",
  weight_id = "weight",
  hvg = NULL,
  scale = TRUE,
  model = c("ns", "std"),
  min_prop = 0.01,
  verbose = TRUE,
  slot_sc = "counts",
  slot_sp = "counts",
  ...
)

Arguments

x, y

single-cell and mixture dataset, respectively. Can be a numeric matrix or a SingleCellExperiment.

groups

vector of group labels for cells in x. When x is a SingleCellExperiment , defaults to colLabels, respectively.

mgs

data.frame or DataFrame of marker genes. Must contain columns holding gene identifiers, group labels and the weight (e.g., logFC, -log(p-value) a feature has in a given group.

n_top

integer scalar specifying the number of markers to select per group. By default NULL uses all the marker genes to initialize the model.

gene_id, group_id, weight_id

character specifying the column in mgs containing gene identifiers, group labels and weights, respectively.

hvg

character vector containing hvg to include in the model. By default NULL.

scale

logical specifying whether to scale single-cell counts to unit variance. This gives the user the option to normalize the data beforehand as you see fit (CPM, FPKM, ...) when passing a matrix or specifying the slot from where to extract the count data.

model

character string indicating which model to use when running NMF. Either "ns" (default) or "std".

min_prop

scalar in [0,1] setting the minimum contribution expected from a cell type in x to observations in y. By default 0.

verbose

logical. Should information on progress be reported?

slot_sc, slot_sp

If the object is of class SingleCellExperiment indicates matrix to use.By default "counts".

...

additional parameters.

Details

SPOTlight uses a Non-Negative Matrix Factorization approach to learn which genes are important for each cell type. In order to drive the factorization and give more importance to cell type marker genes we previously compute them and use them to initialize the basis matrix. This initialized matrices will then be used to carry out the factorization with the single cell expression data. Once the model has learn the topic profiles for each cell type we use non-negative least squares (NNLS) to obtain the topic contributions to each spot. Lastly, NNLS is again used to obtain the proportion of each cell type for each spot by finding the fitting the single-cell topic profiles to the spots topic contributions.

Value

a numeric matrix with rows corresponding to samples and columns to groups

Author(s)

Marc Elosua-Bayes & Helena L. Crowell

Examples

library(scater)
library(scran)

# Use Mock data
# Refer to the vignette for a full workflow
sce <- mockSC(ng = 200, nc = 10, nt = 3)
spe <- mockSP(sce)
mgs <- getMGS(sce)

res <- SPOTlight(
    x = counts(sce),
    y = counts(spe),
    groups = as.character(sce$type),
    mgs = mgs,
    hvg = NULL,
    weight_id = "weight",
    group_id = "type",
    gene_id = "gene")

MarcElosua/SPOTlight documentation built on March 7, 2024, 4:58 p.m.