Chapter 6 Muraro human pancreas (CEL-seq)

6.1 Introduction

This performs an analysis of the Muraro et al. (2016) CEL-seq dataset, consisting of human pancreas cells from various donors.

6.2 Data loading

library(scRNAseq)
sce.muraro <- MuraroPancreasData()

Converting back to Ensembl identifiers.

library(AnnotationHub)
edb <- AnnotationHub()[["AH73881"]]
gene.symb <- sub("__chr.*$", "", rownames(sce.muraro))
gene.ids <- mapIds(edb, keys=gene.symb, 
    keytype="SYMBOL", column="GENEID")

# Removing duplicated genes or genes without Ensembl IDs.
keep <- !is.na(gene.ids) & !duplicated(gene.ids)
sce.muraro <- sce.muraro[keep,]
rownames(sce.muraro) <- gene.ids[keep]

6.3 Quality control

unfiltered <- sce.muraro

This dataset lacks mitochondrial genes so we will do without. For the one batch that seems to have a high proportion of low-quality cells, we compute an appropriate filter threshold using a shared median and MAD from the other batches (Figure 6.1).

library(scater)
stats <- perCellQCMetrics(sce.muraro)
qc <- quickPerCellQC(stats, percent_subsets="altexps_ERCC_percent",
    batch=sce.muraro$donor, subset=sce.muraro$donor!="D28")
sce.muraro <- sce.muraro[,!qc$discard]
colData(unfiltered) <- cbind(colData(unfiltered), stats)
unfiltered$discard <- qc$discard

gridExtra::grid.arrange(
    plotColData(unfiltered, x="donor", y="sum", colour_by="discard") +
        scale_y_log10() + ggtitle("Total count"),
    plotColData(unfiltered, x="donor", y="detected", colour_by="discard") +
        scale_y_log10() + ggtitle("Detected features"),
    plotColData(unfiltered, x="donor", y="altexps_ERCC_percent",
        colour_by="discard") + ggtitle("ERCC percent"),
    ncol=2
)
Distribution of each QC metric across cells from each donor in the Muraro pancreas dataset. Each point represents a cell and is colored according to whether that cell was discarded.

Figure 6.1: Distribution of each QC metric across cells from each donor in the Muraro pancreas dataset. Each point represents a cell and is colored according to whether that cell was discarded.

We have a look at the causes of removal:

colSums(as.matrix(qc))
##              low_lib_size            low_n_features high_altexps_ERCC_percent 
##                       663                       700                       738 
##                   discard 
##                       773

6.4 Normalization

library(scran)
set.seed(1000)
clusters <- quickCluster(sce.muraro)
sce.muraro <- computeSumFactors(sce.muraro, clusters=clusters)
sce.muraro <- logNormCounts(sce.muraro)
summary(sizeFactors(sce.muraro))
##    Min. 1st Qu.  Median    Mean 3rd Qu.    Max. 
##   0.088   0.541   0.821   1.000   1.211  13.987
plot(librarySizeFactors(sce.muraro), sizeFactors(sce.muraro), pch=16,
    xlab="Library size factors", ylab="Deconvolution factors", log="xy")
Relationship between the library size factors and the deconvolution size factors in the Muraro pancreas dataset.

Figure 6.2: Relationship between the library size factors and the deconvolution size factors in the Muraro pancreas dataset.

6.5 Variance modelling

We block on a combined plate and donor factor.

block <- paste0(sce.muraro$plate, "_", sce.muraro$donor)
dec.muraro <- modelGeneVarWithSpikes(sce.muraro, "ERCC", block=block)
top.muraro <- getTopHVGs(dec.muraro, prop=0.1)
par(mfrow=c(8,4))
blocked.stats <- dec.muraro$per.block
for (i in colnames(blocked.stats)) {
    current <- blocked.stats[[i]]
    plot(current$mean, current$total, main=i, pch=16, cex=0.5,
        xlab="Mean of log-expression", ylab="Variance of log-expression")
    curfit <- metadata(current)
    points(curfit$mean, curfit$var, col="red", pch=16)
    curve(curfit$trend(x), col='dodgerblue', add=TRUE, lwd=2)
}
Per-gene variance as a function of the mean for the log-expression values in the Muraro pancreas dataset. Each point represents a gene (black) with the mean-variance trend (blue) fitted to the spike-in transcripts (red) separately for each donor.

Figure 6.3: Per-gene variance as a function of the mean for the log-expression values in the Muraro pancreas dataset. Each point represents a gene (black) with the mean-variance trend (blue) fitted to the spike-in transcripts (red) separately for each donor.

6.6 Data integration

library(batchelor)
set.seed(1001010)
merged.muraro <- fastMNN(sce.muraro, subset.row=top.muraro, 
    batch=sce.muraro$donor)

We use the proportion of variance lost as a diagnostic measure:

metadata(merged.muraro)$merge.info$lost.var
##           D28      D29      D30     D31
## [1,] 0.060847 0.024121 0.000000 0.00000
## [2,] 0.002646 0.003018 0.062421 0.00000
## [3,] 0.003449 0.002641 0.002598 0.08162

6.7 Dimensionality reduction

set.seed(100111)
merged.muraro <- runTSNE(merged.muraro, dimred="corrected")

6.8 Clustering

snn.gr <- buildSNNGraph(merged.muraro, use.dimred="corrected")
colLabels(merged.muraro) <- factor(igraph::cluster_walktrap(snn.gr)$membership)
tab <- table(Cluster=colLabels(merged.muraro), CellType=sce.muraro$label)
library(pheatmap)
pheatmap(log10(tab+10), color=viridis::viridis(100))
Heatmap of the frequency of cells from each cell type label in each cluster.

Figure 6.4: Heatmap of the frequency of cells from each cell type label in each cluster.

table(Cluster=colLabels(merged.muraro), Donor=merged.muraro$batch)
##        Donor
## Cluster D28 D29 D30 D31
##      1  104   6  57 112
##      2   59  21  77  97
##      3   12  75  64  43
##      4   28 149 126 120
##      5   87 261 277 214
##      6   21   7  54  26
##      7    1   6   6  37
##      8    6   6   5   2
##      9   11  68   5  30
##      10   4   2   5   8
gridExtra::grid.arrange(
    plotTSNE(merged.muraro, colour_by="label"),
    plotTSNE(merged.muraro, colour_by="batch"),
    ncol=2
)
Obligatory $t$-SNE plots of the Muraro pancreas dataset. Each point represents a cell that is colored by cluster (left) or batch (right).

Figure 6.5: Obligatory \(t\)-SNE plots of the Muraro pancreas dataset. Each point represents a cell that is colored by cluster (left) or batch (right).

Session Info

R version 4.3.1 (2023-06-16)
Platform: x86_64-pc-linux-gnu (64-bit)
Running under: Ubuntu 22.04.3 LTS

Matrix products: default
BLAS:   /home/biocbuild/bbs-3.18-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_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       

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] pheatmap_1.0.12             BiocSingular_1.18.0        
 [3] batchelor_1.18.0            scran_1.30.0               
 [5] scater_1.30.0               ggplot2_3.4.4              
 [7] scuttle_1.12.0              ensembldb_2.26.0           
 [9] AnnotationFilter_1.26.0     GenomicFeatures_1.54.0     
[11] AnnotationDbi_1.64.0        AnnotationHub_3.10.0       
[13] BiocFileCache_2.10.0        dbplyr_2.3.4               
[15] Matrix_1.6-1.1              scRNAseq_2.15.0            
[17] SingleCellExperiment_1.24.0 SummarizedExperiment_1.32.0
[19] Biobase_2.62.0              GenomicRanges_1.54.0       
[21] GenomeInfoDb_1.38.0         IRanges_2.36.0             
[23] S4Vectors_0.40.0            BiocGenerics_0.48.0        
[25] MatrixGenerics_1.14.0       matrixStats_1.0.0          
[27] BiocStyle_2.30.0            rebook_1.12.0              

loaded via a namespace (and not attached):
  [1] RColorBrewer_1.1-3            rstudioapi_0.15.0            
  [3] jsonlite_1.8.7                CodeDepends_0.6.5            
  [5] magrittr_2.0.3                ggbeeswarm_0.7.2             
  [7] farver_2.1.1                  rmarkdown_2.25               
  [9] BiocIO_1.12.0                 zlibbioc_1.48.0              
 [11] vctrs_0.6.4                   memoise_2.0.1                
 [13] Rsamtools_2.18.0              DelayedMatrixStats_1.24.0    
 [15] RCurl_1.98-1.12               htmltools_0.5.6.1            
 [17] S4Arrays_1.2.0                progress_1.2.2               
 [19] curl_5.1.0                    BiocNeighbors_1.20.0         
 [21] SparseArray_1.2.0             sass_0.4.7                   
 [23] bslib_0.5.1                   cachem_1.0.8                 
 [25] ResidualMatrix_1.12.0         GenomicAlignments_1.38.0     
 [27] igraph_1.5.1                  mime_0.12                    
 [29] lifecycle_1.0.3               pkgconfig_2.0.3              
 [31] rsvd_1.0.5                    R6_2.5.1                     
 [33] fastmap_1.1.1                 GenomeInfoDbData_1.2.11      
 [35] shiny_1.7.5.1                 digest_0.6.33                
 [37] colorspace_2.1-0              dqrng_0.3.1                  
 [39] irlba_2.3.5.1                 ExperimentHub_2.10.0         
 [41] RSQLite_2.3.1                 beachmat_2.18.0              
 [43] labeling_0.4.3                filelock_1.0.2               
 [45] fansi_1.0.5                   httr_1.4.7                   
 [47] abind_1.4-5                   compiler_4.3.1               
 [49] bit64_4.0.5                   withr_2.5.1                  
 [51] BiocParallel_1.36.0           viridis_0.6.4                
 [53] DBI_1.1.3                     biomaRt_2.58.0               
 [55] rappdirs_0.3.3                DelayedArray_0.28.0          
 [57] bluster_1.12.0                rjson_0.2.21                 
 [59] tools_4.3.1                   vipor_0.4.5                  
 [61] beeswarm_0.4.0                interactiveDisplayBase_1.40.0
 [63] httpuv_1.6.12                 glue_1.6.2                   
 [65] restfulr_0.0.15               promises_1.2.1               
 [67] grid_4.3.1                    Rtsne_0.16                   
 [69] cluster_2.1.4                 generics_0.1.3               
 [71] gtable_0.3.4                  hms_1.1.3                    
 [73] metapod_1.10.0                ScaledMatrix_1.10.0          
 [75] xml2_1.3.5                    utf8_1.2.4                   
 [77] XVector_0.42.0                ggrepel_0.9.4                
 [79] BiocVersion_3.18.0            pillar_1.9.0                 
 [81] stringr_1.5.0                 limma_3.58.0                 
 [83] later_1.3.1                   dplyr_1.1.3                  
 [85] lattice_0.22-5                rtracklayer_1.62.0           
 [87] bit_4.0.5                     tidyselect_1.2.0             
 [89] locfit_1.5-9.8                Biostrings_2.70.0            
 [91] knitr_1.44                    gridExtra_2.3                
 [93] bookdown_0.36                 ProtGenerics_1.34.0          
 [95] edgeR_4.0.0                   xfun_0.40                    
 [97] statmod_1.5.0                 stringi_1.7.12               
 [99] lazyeval_0.2.2                yaml_2.3.7                   
[101] evaluate_0.22                 codetools_0.2-19             
[103] tibble_3.2.1                  BiocManager_1.30.22          
[105] graph_1.80.0                  cli_3.6.1                    
[107] xtable_1.8-4                  munsell_0.5.0                
[109] jquerylib_0.1.4               Rcpp_1.0.11                  
[111] dir.expiry_1.10.0             png_0.1-8                    
[113] XML_3.99-0.14                 parallel_4.3.1               
[115] ellipsis_0.3.2                blob_1.2.4                   
[117] prettyunits_1.2.0             sparseMatrixStats_1.14.0     
[119] bitops_1.0-7                  viridisLite_0.4.2            
[121] scales_1.2.1                  purrr_1.0.2                  
[123] crayon_1.5.2                  rlang_1.1.1                  
[125] cowplot_1.1.1                 KEGGREST_1.42.0              

References

Muraro, M. J., G. Dharmadhikari, D. Grun, N. Groen, T. Dielen, E. Jansen, L. van Gurp, et al. 2016. “A Single-Cell Transcriptome Atlas of the Human Pancreas.” Cell Syst 3 (4): 385–94.