scTGIFscTGIF 1.18.0
Here, we explain the concept of scTGIF. The analysis of single-cell RNA-Seq (scRNA-Seq) has a potential difficult problem; which data corresponds to what kind of cell type is not known a priori.
Therefore, at the start point of the data analysis of the scRNA-Seq dataset, each cell is “not colored” (unannotated) (Figure 1). There some approaches to support users to infer the cell types such as (1) Known marker gene expression, (2) BLAST-like gene expression comparison with reference DB, (3) differentially expressed genes (DEGs) and over-representative analysis (ORA) (scRNA-tools).
The first approach might be the most popular method, but this task is based on the expert knowledge about the cell types, and not always general-purpose. The second approach is easy and scalable, but still limited when the cell type is not known or still not measured by the other research organization. The third approach can perhaps be used in any situation but ambiguous and time-consuming task; this task is based on the cluster label and the true cluster structure, which is not known and some DEG methods have to be performed in each cluster, but recent scRNA-Seq dataset has tens to hundreds of cell types. Besides, a scRNA-Seq dataset can have low-quality cells and artifacts (e.g. doublet) but it is hard to distinguish from real cell data. Therefore, in actual data analytical situation, laborious trial-and-error cycle along with the change of cellular label cannot be evitable (Figure 1).
scTGIF is developed to reduce this trial-and-error cycle; This tool directly connects the unannotated cells and related gene function. Since this tool does not use reference DB, marker gene list, and cluster label can be used in any situation without expert knowledge and is not influenced by the change of cellular label.
Figure 1: Concept of scTGIF
In scTGIF, three data is required; the gene expression matrix, 2D coordinates of the cells (e.g. t-SNE, UMAP), and geneset of MSigDB. Firstly, the 2D coordinates are segmented as 50-by-50 grids, and gene expression is summarized in each grid level (X1). Next, the correspondence between genes and the related gene functions are summarized as gene-by-function matrix (X2). Here, we support only common genes are used in X1 and X2. Performing joint non-negative matrix factorization (jNMF) algorithm, which is implemented in nnTensor, the shared latent variables (W) with the two matrices are estimated.
Figure 2: Joint NMF
By this algorithm, a grid set and corresponding gene functions are paired. Lower-dimension (D)-by-Grid matrix H1 works as attention maps to help users to pay attention the grids, and D-by-Function matrix H2 shows the gene function enriched in the grids.
Figure 3: H1 and H2 matrices
scTGIF also supports some QC metrics to distinguish low-quality cells and artifacts from real cellular data.
To demonstrate the usage of scTGIF, we prepared a testdata of distal lung epithelium.
library("scTGIF")
library("SingleCellExperiment")
library("GSEABase")
library("msigdbr")
data("DistalLungEpithelium")
data("pca.DistalLungEpithelium")
data("label.DistalLungEpithelium")Although this data is still annotated and the cell type label is provided, scTGIF does not rely on this information.
par(ask=FALSE)
plot(pca.DistalLungEpithelium, col=label.DistalLungEpithelium, pch=16,
    main="Distal lung epithelium dataset", xlab="PCA1", ylab="PCA2", bty="n")
text(0.1, 0.05, "AT1", col="#FF7F00", cex=2)
text(0.07, -0.15, "AT2", col="#E41A1C", cex=2)
text(0.13, -0.04, "BP", col="#A65628", cex=2)
text(0.125, -0.15, "Clara", col="#377EB8", cex=2)
text(0.09, -0.2, "Cilliated", col="#4DAF4A", cex=2)To combine with the gene expression and the related gene function, we suppose the gene function data is summarized as the object of GSEABase. This data is directly downloadable from MSigDB and can be imported like gmt <- GSEABase::getGmt(“/YOURPATH/h.all.v6.0.entrez.gmt”)
Note that scTGIF only supports NCBI Gene IDs (Entrez IDs) for now. When the scRNA-Seq is not about human, the situation is more complicated. Here, we use msigdbr package to retrieve the mouse MSigDB genesets like below.
m_df = msigdbr(species = "Mus musculus",
    category = "H")[, c("gs_name", "entrez_gene")]
hallmark = unique(m_df$gs_name)
gsc <- lapply(hallmark, function(h){
    target = which(m_df$gs_name == h)
    geneIds = unique(as.character(m_df$entrez_gene[target]))
    GeneSet(setName=h, geneIds)
})
gmt = GeneSetCollection(gsc)
gmt = gmt[1:10] # Reduced for this demoNext, the data matrix is converted to the object of SingleCellExperiment package, and the 2D coordinates are registered as the reducedDims slot.
sce <- SingleCellExperiment(assays = list(counts = DistalLungEpithelium))
reducedDims(sce) <- SimpleList(PCA=pca.DistalLungEpithelium)Although the default mode of scTGIF use count slot as
the input matrix,
the normalized gene expression data can also be specified.
In such a case, we recommend using normcounts slot to register the data.
CPMED <- function(input){
    libsize <- colSums(input)
    median(libsize) * t(t(input) / libsize)
}
normcounts(sce) <- log10(CPMED(counts(sce)) + 1)After the registration of the data in sce,
settingTGIF will work like below.
settingTGIF(sce, gmt, reducedDimNames="PCA", assayNames="normcounts")Finally, reportTGIF generates the HTML report to summarize
the result of jNMF.
reportTGIF(sce,
    html.open=FALSE,
    title="scTGIF Report for DistalLungEpithelium dataset",
    author="Koki Tsuyuzaki")## index.Rmd is created...## index.Rmd is compiled to index.html...## ################################################
## Data files are saved in
## /tmp/Rtmpc86CxE
## ################################################Since this function takes some time,
please type example("reportTGIF") by your own environment.
## R version 4.4.0 beta (2024-04-15 r86425)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 22.04.4 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.19-bioc/R/lib/libRblas.so 
## LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.10.0
## 
## locale:
##  [1] LC_CTYPE=en_US.UTF-8       LC_NUMERIC=C              
##  [3] LC_TIME=en_US.UTF-8        LC_COLLATE=en_US.UTF-8    
##  [5] LC_MONETARY=en_US.UTF-8    LC_MESSAGES=en_US.UTF-8   
##  [7] LC_PAPER=en_US.UTF-8       LC_NAME=C                 
##  [9] LC_ADDRESS=C               LC_TELEPHONE=C            
## [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C       
## 
## time zone: America/New_York
## tzcode source: system (glibc)
## 
## attached base packages:
## [1] stats4    stats     graphics  grDevices utils     datasets  methods  
## [8] base     
## 
## other attached packages:
##  [1] plotly_4.10.4               ggplot2_3.5.1              
##  [3] msigdbr_7.5.1               GSEABase_1.66.0            
##  [5] graph_1.82.0                annotate_1.82.0            
##  [7] XML_3.99-0.16.1             AnnotationDbi_1.66.0       
##  [9] SingleCellExperiment_1.26.0 SummarizedExperiment_1.34.0
## [11] Biobase_2.64.0              GenomicRanges_1.56.0       
## [13] GenomeInfoDb_1.40.0         IRanges_2.38.0             
## [15] S4Vectors_0.42.0            BiocGenerics_0.50.0        
## [17] MatrixGenerics_1.16.0       matrixStats_1.3.0          
## [19] scTGIF_1.18.0               BiocStyle_2.32.0           
## 
## loaded via a namespace (and not attached):
##  [1] DBI_1.2.2               tcltk_4.4.0             rlang_1.1.3            
##  [4] magrittr_2.0.3          compiler_4.4.0          RSQLite_2.3.6          
##  [7] systemfonts_1.0.6       png_0.1-8               vctrs_0.6.5            
## [10] maps_3.4.2              nnTensor_1.2.0          pkgconfig_2.0.3        
## [13] crayon_1.5.2            fastmap_1.1.1           magick_2.8.3           
## [16] XVector_0.44.0          labeling_0.4.3          utf8_1.2.4             
## [19] tagcloud_0.6            rmarkdown_2.26          UCSC.utils_1.0.0       
## [22] ragg_1.3.0              tinytex_0.50            purrr_1.0.2            
## [25] bit_4.0.5               xfun_0.43               zlibbioc_1.50.0        
## [28] cachem_1.0.8            jsonlite_1.8.8          blob_1.2.4             
## [31] highr_0.10              DelayedArray_0.30.0     tweenr_2.0.3           
## [34] cluster_2.1.6           R6_2.5.1                bslib_0.7.0            
## [37] RColorBrewer_1.1-3      schex_1.18.0            jquerylib_0.1.4        
## [40] Rcpp_1.0.12             bookdown_0.39           knitr_1.46             
## [43] fields_15.2             igraph_2.0.3            Matrix_1.7-0           
## [46] tidyselect_1.2.1        abind_1.4-5             yaml_2.3.8             
## [49] misc3d_0.9-1            lattice_0.22-6          tibble_3.2.1           
## [52] withr_3.0.0             KEGGREST_1.44.0         evaluate_0.23          
## [55] polyclip_1.10-6         Biostrings_2.72.0       pillar_1.9.0           
## [58] BiocManager_1.30.22     generics_0.1.3          munsell_0.5.1          
## [61] scales_1.3.0            xtable_1.8-4            rTensor_1.4.8          
## [64] glue_1.7.0              lazyeval_0.2.2          tools_4.4.0            
## [67] hexbin_1.28.3           data.table_1.15.4       babelgene_22.9         
## [70] dotCall64_1.1-1         grid_4.4.0              tidyr_1.3.1            
## [73] crosstalk_1.2.1         colorspace_2.1-0        GenomeInfoDbData_1.2.12
## [76] ggforce_0.4.2           cli_3.6.2               textshaping_0.3.7      
## [79] spam_2.10-0             fansi_1.0.6             S4Arrays_1.4.0         
## [82] plot3D_1.4.1            viridisLite_0.4.2       concaveman_1.1.0       
## [85] dplyr_1.1.4             gtable_0.3.5            sass_0.4.9             
## [88] digest_0.6.35           SparseArray_1.4.0       farver_2.1.1           
## [91] htmlwidgets_1.6.4       entropy_1.3.1           memoise_2.0.1          
## [94] htmltools_0.5.8.1       lifecycle_1.0.4         httr_1.4.7             
## [97] bit64_4.0.5             MASS_7.3-60.2