1 Introduction

snifter provides an R wrapper for the openTSNE implementation of fast interpolated t-SNE (FI-tSNE). It is based on basilisk and reticulate. This vignette aims to provide a brief overview of typical use when applied to scRNAseq data, but it does not provide a comprehensive guide to the available options in the package.

It is highly advisable to review the documentation in snifter and the openTSNE documentation to gain a full understanding of the available options.

2 Setting up the data

We will illustrate the use of snifter using data from scRNAseq and single cell utility functions provided by scuttle, scater and scran - first we load these libraries and set a random seed to ensure the t-SNE visualisation is reproducible (note: it is good practice to ensure that a t-SNE embedding is robust by running the algorithm multiple times).

library("snifter")
library("scRNAseq")
library("scran")
library("scuttle")
library("scater")
library("ggplot2")
theme_set(theme_bw())
set.seed(42)

Before running t-SNE, we first load data generated by Zeisel et al. from scRNAseq. We filter this data to remove genes expressed only in a small number of cells, estimate normalisation factors using scran and generate 20 principal components. We will use these principal components to generate the t-SNE embedding later.

data <- ZeiselBrainData()
data <- data[rowMeans(counts(data) != 0) > 0.05, ]
data <- computeSumFactors(data, cluster = quickCluster(data))
data <- logNormCounts(data)
data <- runPCA(data, ncomponents = 20)
## Convert this to a factor to use as colouring variable later
data$level1class <- factor(data$level1class)

3 Running t-SNE

The main functionality of the package lies in the fitsne function. This function returns a matrix of t-SNE co-ordinates. In this case, we pass in the 20 principal components computed based on the log-normalised counts. We colour points based on the discrete cell types identified by the authors.

mat <- reducedDim(data)
fit <- fitsne(mat, random_state = 42L)
ggplot() +
    aes(fit[, 1], fit[, 2], colour = data$level1class) +
    geom_point(pch = 19) +
    scale_colour_discrete(name = "Cell type") +
    labs(x = "t-SNE 1", y = "t-SNE 2")

4 Projecting new data into an existing embedding

The openTNSE package, and by extension snifter, also allows the embedding of new data into an existing t-SNE embedding. Here, we will split the data into “training” and “test” sets. Following this, we generate a t-SNE embedding using the training data, and project the test data into this embedding.

test_ind <- sample(nrow(mat), nrow(mat) / 2)
train_ind <- setdiff(seq_len(nrow(mat)), test_ind)
train_mat <- mat[train_ind, ]
test_mat <- mat[test_ind, ]

train_label <- data$level1class[train_ind]
test_label <- data$level1class[test_ind]

embedding <- fitsne(train_mat, random_state = 42L)

Once we have generated the embedding, we can now project the unseen test data into this t-SNE embedding.

new_coords <- project(embedding, new = test_mat, old = train_mat)
ggplot() +
    geom_point(
        aes(embedding[, 1], embedding[, 2],
            colour = train_label,
            shape = "Train"
        )
    ) +
    geom_point(
        aes(new_coords[, 1], new_coords[, 2], 
            colour = test_label,
            shape = "Test"
        )
    ) +
    scale_colour_discrete(name = "Cell type") +
    scale_shape_discrete(name = NULL) +
    labs(x = "t-SNE 1", y = "t-SNE 2")

Session information

sessionInfo()
#> R version 4.1.1 (2021-08-10)
#> Platform: x86_64-pc-linux-gnu (64-bit)
#> Running under: Ubuntu 20.04.3 LTS
#> 
#> Matrix products: default
#> BLAS:   /home/biocbuild/bbs-3.14-bioc/R/lib/libRblas.so
#> LAPACK: /home/biocbuild/bbs-3.14-bioc/R/lib/libRlapack.so
#> 
#> locale:
#>  [1] LC_CTYPE=en_US.UTF-8       LC_NUMERIC=C              
#>  [3] LC_TIME=en_GB              LC_COLLATE=C              
#>  [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       
#> 
#> attached base packages:
#> [1] stats4    stats     graphics  grDevices utils     datasets  methods  
#> [8] base     
#> 
#> other attached packages:
#>  [1] scater_1.22.0               ggplot2_3.3.5              
#>  [3] scran_1.22.0                scuttle_1.4.0              
#>  [5] scRNAseq_2.7.2              SingleCellExperiment_1.16.0
#>  [7] SummarizedExperiment_1.24.0 Biobase_2.54.0             
#>  [9] GenomicRanges_1.46.0        GenomeInfoDb_1.30.0        
#> [11] IRanges_2.28.0              S4Vectors_0.32.0           
#> [13] BiocGenerics_0.40.0         MatrixGenerics_1.6.0       
#> [15] matrixStats_0.61.0          snifter_1.4.0              
#> [17] BiocStyle_2.22.0           
#> 
#> loaded via a namespace (and not attached):
#>   [1] AnnotationHub_3.2.0           BiocFileCache_2.2.0          
#>   [3] igraph_1.2.7                  lazyeval_0.2.2               
#>   [5] BiocParallel_1.28.0           digest_0.6.28                
#>   [7] ensembldb_2.18.0              htmltools_0.5.2              
#>   [9] magick_2.7.3                  viridis_0.6.2                
#>  [11] fansi_0.5.0                   magrittr_2.0.1               
#>  [13] memoise_2.0.0                 ScaledMatrix_1.2.0           
#>  [15] cluster_2.1.2                 limma_3.50.0                 
#>  [17] Biostrings_2.62.0             prettyunits_1.1.1            
#>  [19] colorspace_2.0-2              ggrepel_0.9.1                
#>  [21] blob_1.2.2                    rappdirs_0.3.3               
#>  [23] xfun_0.27                     dplyr_1.0.7                  
#>  [25] crayon_1.4.1                  RCurl_1.98-1.5               
#>  [27] jsonlite_1.7.2                glue_1.4.2                   
#>  [29] gtable_0.3.0                  zlibbioc_1.40.0              
#>  [31] XVector_0.34.0                DelayedArray_0.20.0          
#>  [33] BiocSingular_1.10.0           scales_1.1.1                 
#>  [35] DBI_1.1.1                     edgeR_3.36.0                 
#>  [37] Rcpp_1.0.7                    viridisLite_0.4.0            
#>  [39] xtable_1.8-4                  progress_1.2.2               
#>  [41] reticulate_1.22               dqrng_0.3.0                  
#>  [43] bit_4.0.4                     rsvd_1.0.5                   
#>  [45] metapod_1.2.0                 httr_1.4.2                   
#>  [47] dir.expiry_1.2.0              ellipsis_0.3.2               
#>  [49] farver_2.1.0                  pkgconfig_2.0.3              
#>  [51] XML_3.99-0.8                  sass_0.4.0                   
#>  [53] dbplyr_2.1.1                  here_1.0.1                   
#>  [55] locfit_1.5-9.4                utf8_1.2.2                   
#>  [57] labeling_0.4.2                tidyselect_1.1.1             
#>  [59] rlang_0.4.12                  later_1.3.0                  
#>  [61] AnnotationDbi_1.56.0          munsell_0.5.0                
#>  [63] BiocVersion_3.14.0            tools_4.1.1                  
#>  [65] cachem_1.0.6                  generics_0.1.1               
#>  [67] RSQLite_2.2.8                 ExperimentHub_2.2.0          
#>  [69] evaluate_0.14                 stringr_1.4.0                
#>  [71] fastmap_1.1.0                 yaml_2.2.1                   
#>  [73] knitr_1.36                    bit64_4.0.5                  
#>  [75] purrr_0.3.4                   KEGGREST_1.34.0              
#>  [77] AnnotationFilter_1.18.0       sparseMatrixStats_1.6.0      
#>  [79] mime_0.12                     xml2_1.3.2                   
#>  [81] biomaRt_2.50.0                compiler_4.1.1               
#>  [83] beeswarm_0.4.0                filelock_1.0.2               
#>  [85] curl_4.3.2                    png_0.1-7                    
#>  [87] interactiveDisplayBase_1.32.0 tibble_3.1.5                 
#>  [89] statmod_1.4.36                bslib_0.3.1                  
#>  [91] stringi_1.7.5                 highr_0.9                    
#>  [93] basilisk.utils_1.6.0          GenomicFeatures_1.46.0       
#>  [95] lattice_0.20-45               bluster_1.4.0                
#>  [97] ProtGenerics_1.26.0           Matrix_1.3-4                 
#>  [99] vctrs_0.3.8                   pillar_1.6.4                 
#> [101] lifecycle_1.0.1               BiocManager_1.30.16          
#> [103] jquerylib_0.1.4               BiocNeighbors_1.12.0         
#> [105] bitops_1.0-7                  irlba_2.3.3                  
#> [107] httpuv_1.6.3                  rtracklayer_1.54.0           
#> [109] R6_2.5.1                      BiocIO_1.4.0                 
#> [111] bookdown_0.24                 promises_1.2.0.1             
#> [113] gridExtra_2.3                 vipor_0.4.5                  
#> [115] assertthat_0.2.1              rprojroot_2.0.2              
#> [117] rjson_0.2.20                  withr_2.4.2                  
#> [119] GenomicAlignments_1.30.0      Rsamtools_2.10.0             
#> [121] GenomeInfoDbData_1.2.7        parallel_4.1.1               
#> [123] hms_1.1.1                     grid_4.1.1                   
#> [125] beachmat_2.10.0               basilisk_1.6.0               
#> [127] rmarkdown_2.11                DelayedMatrixStats_1.16.0    
#> [129] shiny_1.7.1                   ggbeeswarm_0.6.0             
#> [131] restfulr_0.0.13