SmoothKNN.Seurat | R Documentation |
This function performs smoothing of single-cell scores by weighted average of the k-nearest neighbors. It can be useful to 'impute' scores by neighboring cells and partially correct data sparsity. While this function has been designed to smooth UCell scores, it can be applied to any numerical metadata contained in SingleCellExperiment or Seurat objects
## S3 method for class 'Seurat'
SmoothKNN(
obj = NULL,
signature.names = NULL,
reduction = "pca",
k = 10,
decay = 0.1,
up.only = FALSE,
BNPARAM = AnnoyParam(),
BPPARAM = SerialParam(),
suffix = "_kNN",
assay = NULL,
slot = "data",
sce.expname = NULL,
sce.assay = NULL
)
## S3 method for class 'SingleCellExperiment'
SmoothKNN(
obj = NULL,
signature.names = NULL,
reduction = "PCA",
k = 10,
decay = 0.1,
up.only = FALSE,
BNPARAM = AnnoyParam(),
BPPARAM = SerialParam(),
suffix = "_kNN",
assay = NULL,
slot = "data",
sce.expname = c("UCell", "main"),
sce.assay = NULL
)
SmoothKNN(
obj = NULL,
signature.names = NULL,
reduction = "pca",
k = 10,
decay = 0.1,
up.only = FALSE,
BNPARAM = AnnoyParam(),
BPPARAM = SerialParam(),
suffix = "_kNN",
assay = NULL,
slot = "data",
sce.expname = c("UCell", "main"),
sce.assay = NULL
)
obj |
Input object - either a SingleCellExperiment object or a Seurat object. |
signature.names |
The names of the signatures (or any numeric metadata column) for which to calculate kNN-smoothed scores |
reduction |
Which dimensionality reduction to use for kNN smoothing. It must be already present in the input object. |
k |
Number of neighbors for kNN smoothing |
decay |
Exponential decay for nearest neighbor weight: (1-decay)^n |
up.only |
If set to TRUE, smoothed scores will only be allowed to increase by smoothing |
BNPARAM |
A BiocNeighborParam object specifying the algorithm to use for kNN calculation. |
BPPARAM |
A |
suffix |
Suffix to append to metadata columns for the new knn-smoothed scores |
assay |
For Seurat objects only - do smoothing on expression data from this assay. When NULL, only looks in metadata |
slot |
For Seurat objects only - do smoothing on expression data from this slot |
sce.expname |
For sce objects only - which experiment stores the signatures to be smoothed. Set to 'main' for smoothing gene expression stored in the main sce experiment. |
sce.assay |
For sce objects only - pull data from this assay |
An augmented obj
with the smoothed signatures. If obj
is a Seurat object, smoothed signatures are added to metadata; if
obj
is a SingleCellExperiment object, smoothed signatures are
returned in a new altExp. See the examples below.
#### Using Seurat ####
library(Seurat)
gene.sets <- list(Tcell = c("CD2","CD3E","CD3D"),
Myeloid = c("SPI1","FCER1G","CSF1R"))
data(sample.matrix)
obj <- Seurat::CreateSeuratObject(sample.matrix)
# Calculate UCell scores
obj <- AddModuleScore_UCell(obj,features = gene.sets, name=NULL)
# Run PCA
obj <- FindVariableFeatures(obj) |> NormalizeData() |> ScaleData() |> RunPCA(npcs=5)
# Smooth signatures
obj <- SmoothKNN(obj, k=3, signature.names=names(gene.sets))
head(obj[[]])
#### Using SingleCellExperiment ####
library(SingleCellExperiment)
library(scater)
data(sample.matrix)
sce <- SingleCellExperiment(list(counts=sample.matrix))
gene.sets <- list( Tcell = c("CD2","CD3E","CD3D"),
Myeloid = c("SPI1","FCER1G","CSF1R"))
# Calculate UCell scores
sce <- ScoreSignatures_UCell(sce, features=gene.sets, name=NULL)
# Run PCA
sce <- logNormCounts(sce)
sce <- runPCA(sce, scale=TRUE, ncomponents=5)
# Smooth signatures
sce <- SmoothKNN(sce, k=3, signature.names=names(gene.sets))
# See results
altExp(sce, 'UCell')
assays(altExp(sce, 'UCell'))
# Plot on UMAP
sce <- runUMAP(sce, dimred="PCA")
plotUMAP(sce, colour_by = "Tcell_kNN", by_exprs_values = "UCell_kNN")
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.