--- title: "Using UCell with Seurat" author: - name: Massimo Andreatta affiliation: Ludwig Institute for Cancer Research, Lausanne Branch, and Department of Oncology, CHUV and University of Lausanne, Epalinges 1066, Switzerland; and Swiss Institute of Bioinformatics, Lausanne, Switzerland - name: Santiago J. Carmona affiliation: Ludwig Institute for Cancer Research, Lausanne Branch, and Department of Oncology, CHUV and University of Lausanne, Epalinges 1066, Switzerland; and Swiss Institute of Bioinformatics, Lausanne, Switzerland output: BiocStyle::html_document: toc_float: true BiocStyle::pdf_document: default package: UCell vignette: | %\VignetteIndexEntry{3. Using UCell with Seurat} %\VignetteEngine{knitr::rmarkdown} %\VignetteEncoding{UTF-8} --- # Introduction The function `AddModuleScore_UCell()` allows operating directly on Seurat objects. UCell scores are calculated from raw counts or normalized data, and returned as metadata columns. The example below defines some simple signatures, and applies them on single-cell data stored in a Seurat object. To see how this function differs from Seurat's own `AddModuleScore()` (not based on per-cell ranks) see [this vignette](https://carmonalab.github.io/UCell_demo/UCell_Seurat_vignette.html). # Get some testing data For this demo, we will download a single-cell dataset of lung cancer ([Zilionis et al. (2019) Immunity](https://pubmed.ncbi.nlm.nih.gov/30979687/)) through the [scRNA-seq](https://bioconductor.org/packages/3.15/data/experiment/html/scRNAseq.html) package. This dataset contains >170,000 single cells; for the sake of simplicity, in this demo will we focus on immune cells, according to the annotations by the authors, and downsample to 5000 cells. ```{r message=F, warning=F, results=F} library(scRNAseq) lung <- ZilionisLungData() immune <- lung$Used & lung$used_in_NSCLC_immune lung <- lung[,immune] lung <- lung[,1:5000] exp.mat <- Matrix::Matrix(counts(lung),sparse = TRUE) ``` # Define gene signatures Here we define some simple gene sets based on the "Human Cell Landscape" signatures [Han et al. (2020) Nature](https://www.nature.com/articles/s41586-020-2157-4). You may edit existing signatures, or add new one as elements in a list. ```{r} signatures <- list( Tcell = c("CD3D","CD3E","CD3G","CD2","TRAC"), Myeloid = c("CD14","LYZ","CSF1R","FCER1G","SPI1","LCK-"), NK = c("KLRD1","NCAM1","NKG7","CD3D-","CD3E-"), Plasma_cell = c("IGKC","IGHG3","IGHG1","IGHA1","CD19-") ) ``` # Run UCell on Seurat object ```{r message=F, warning=F} library(UCell) library(Seurat) seurat.object <- CreateSeuratObject(counts = exp.mat, project = "Zilionis_immune") seurat.object <- AddModuleScore_UCell(seurat.object, features=signatures, name=NULL) head(seurat.object[[]]) ``` Generate PCA and UMAP embeddings ```{r message=F, warning=F} seurat.object <- NormalizeData(seurat.object) seurat.object <- FindVariableFeatures(seurat.object, selection.method = "vst", nfeatures = 500) seurat.object <- ScaleData(seurat.object) seurat.object <- RunPCA(seurat.object, npcs = 20, features=VariableFeatures(seurat.object)) seurat.object <- RunUMAP(seurat.object, reduction = "pca", dims = 1:20, seed.use=123) ``` Visualize UCell scores on low-dimensional representation (UMAP) ```{r fig.width=9, fig.height=7, dpi=60} library(ggplot2) library(patchwork) FeaturePlot(seurat.object, reduction = "umap", features = names(signatures)) ``` # Signature smoothing Single-cell data are sparse. It can be useful to 'impute' scores by neighboring cells and partially correct this sparsity. The function `SmoothKNN` performs smoothing of single-cell scores by weighted average of the k-nearest neighbors in a given dimensionality reduction. It can be applied directly on Seurat objects to smooth UCell scores: ```{r} seurat.object <- SmoothKNN(seurat.object, signature.names = names(signatures), reduction="pca") ``` ```{r fig.wide=TRUE, dpi=60} FeaturePlot(seurat.object, reduction = "umap", features = c("NK","NK_kNN")) ``` Smoothing (or imputation) has been designed for UCell scores, but it can be applied to any other data or metadata. For instance, we can perform knn-smoothing directly on gene expression measurements: ```{r warning=FALSE, fig.width=9, fig.height=7, dpi=60} genes <- c("CD2","CSF1R") seurat.object <- SmoothKNN(seurat.object, signature.names=genes, assay="RNA", reduction="pca", k=20, suffix = "_smooth") DefaultAssay(seurat.object) <- "RNA" a <- FeaturePlot(seurat.object, reduction = "umap", features = genes) DefaultAssay(seurat.object) <- "RNA_smooth" b <- FeaturePlot(seurat.object, reduction = "umap", features = genes) a / b ``` # Resources Please report any issues at the [UCell GitHub repository](https://github.com/carmonalab/UCell). More demos available on the [Bioc landing page](https://bioconductor.org/packages/release/bioc/html/UCell.html) and at the [UCell demo repository](https://github.com/carmonalab/UCell_demo). If you find UCell useful, you may also check out the [scGate package](https://github.com/carmonalab/scGate), which relies on UCell scores to automatically purify populations of interest based on gene signatures. See also [SignatuR](https://github.com/carmonalab/SignatuR) for easy storing and retrieval of gene signatures. # References * Andreatta, M., Carmona, S. J. (2021) *UCell: Robust and scalable single-cell gene signature scoring* Computational and Structural Biotechnology Journal * Zilionis, R., Engblom, C., ..., Klein, A. M. (2019) *Single-Cell Transcriptomics of Human and Mouse Lung Cancers Reveals Conserved Myeloid Populations across Individuals and Species* Immunity * Hao, Yuhan, et al. (2021) *Integrated analysis of multimodal single-cell data* Cell # Session Info ```{r} sessionInfo() ```