spatialdecon: Mixed cell deconvolution of spatiall-resolved gene expression...

View source: R/spatialdecon.R

spatialdeconR Documentation

Mixed cell deconvolution of spatiall-resolved gene expression data

Description

Runs the spatialdecon algorithm with added optional functionalities. Workflow is:

  1. compute weights from raw data

  2. Estimate a tumor profile and merge it into the cell profiles matrix

  3. run deconvolution once

  4. remove poorly-fit genes from first round of decon

  5. re-run decon with cleaned-up gene set

  6. combine closely-related cell types

  7. compute p-values

  8. rescale abundance estimates, to proportions of total, proportions of immune, cell counts

Usage

spatialdecon(
  norm,
  bg,
  X = NULL,
  raw = NULL,
  wts = NULL,
  resid_thresh = 3,
  lower_thresh = 0.5,
  align_genes = TRUE,
  is_pure_tumor = NULL,
  n_tumor_clusters = 10,
  cell_counts = NULL,
  cellmerges = NULL,
  maxit = 1000
)

Arguments

norm

p-length expression vector or p * N expression matrix - the actual (linear-scale) data

bg

Same dimension as norm: the background expected at each data point.

X

Cell profile matrix. If NULL, the safeTME matrix is used.

raw

Optional for using an error model to weight the data points. p-length expression vector or p * N expression matrix - the raw (linear-scale) data

wts

Optional, a matrix of weights.

resid_thresh

A scalar, sets a threshold on how extreme individual data points' values can be (in log2 units) before getting flagged as outliers and set to NA.

lower_thresh

A scalar. Before log2-scale residuals are calculated, both observed and fitted values get thresholded up to this value. Prevents log2-scale residuals from becoming extreme in points near zero.

align_genes

Logical. If TRUE, then Y, X, bg, and wts are row-aligned by shared genes.

is_pure_tumor

A logical vector denoting whether each AOI consists of pure tumor. If specified, then the algorithm will derive a tumor expression profile and merge it with the immune profiles matrix.

n_tumor_clusters

Number of tumor-specific columns to merge into the cell profile matrix. Has an impact only when is_pure_tumor argument is used to indicate pure tumor AOIs. Takes this many clusters from the pure-tumor AOI data and gets the average expression profile in each cluster. Default 10.

cell_counts

Number of cells estimated to be within each sample. If provided alongside norm_factors, then the algorithm will additionally output cell abundance esimtates on the scale of cell counts.

cellmerges

A list object holding the mapping from beta's cell names to combined cell names. If left NULL, then defaults to a mapping of granular immune cell definitions to broader categories.

maxit

Maximum number of iterations. Default 1000.

Value

a list:

  • beta: matrix of cell abundance estimates, cells in rows and observations in columns

  • sigmas: covariance matrices of each observation's beta estimates

  • p: matrix of p-values for H0: beta == 0

  • t: matrix of t-statistics for H0: beta == 0

  • se: matrix of standard errors of beta values

  • prop_of_all: rescaling of beta to sum to 1 in each observation

  • prop_of_nontumor: rescaling of beta to sum to 1 in each observation, excluding tumor abundance estimates

  • cell.counts: beta rescaled to estimate cell numbers, based on prop_of_all and nuclei count

  • beta.granular: cell abundances prior to combining closely-related cell types

  • sigma.granular: sigmas prior to combining closely-related cell types

  • cell.counts.granular: cell.counts prior to combining closely-related cell types

  • resids: a matrix of residuals from the model fit. (log2(pmax(y, lower_thresh)) - log2(pmax(xb, lower_thresh))).

  • X: the cell profile matrix used in the decon fit.

Examples

data(mini_geomx_dataset)
data(safeTME)
data(safeTME.matches)
# estimate background:
mini_geomx_dataset$bg <- derive_GeoMx_background(
  norm = mini_geomx_dataset$normalized,
  probepool = rep(1, nrow(mini_geomx_dataset$normalized)),
  negnames = "NegProbe"
)
# run basic decon:
res0 <- spatialdecon(
  norm = mini_geomx_dataset$normalized,
  bg = mini_geomx_dataset$bg,
  X = safeTME
)
# run decon with bells and whistles:
res <- spatialdecon(
  norm = mini_geomx_dataset$normalized,
  bg = mini_geomx_dataset$bg,
  X = safeTME,
  cellmerges = safeTME.matches,
  cell_counts = mini_geomx_dataset$annot$nuclei,
  is_pure_tumor = mini_geomx_dataset$annot$AOI.name == "Tumor"
)

Nanostring-Biostats/SpatialDecon documentation built on Jan. 26, 2024, 8:20 p.m.