Chapter 5 Grun human pancreas (CEL-seq2)

5.1 Introduction

This workflow performs an analysis of the Grun et al. (2016) CEL-seq2 dataset consisting of human pancreas cells from various donors.

5.2 Data loading

library(scRNAseq)
sce.grun <- GrunPancreasData()

We convert to Ensembl identifiers, and we remove duplicated genes or genes without Ensembl IDs.

library(org.Hs.eg.db)
gene.ids <- mapIds(org.Hs.eg.db, keys=rowData(sce.grun)$symbol,
    keytype="SYMBOL", column="ENSEMBL")

keep <- !is.na(gene.ids) & !duplicated(gene.ids)
sce.grun <- sce.grun[keep,]
rownames(sce.grun) <- gene.ids[keep]

5.3 Quality control

unfiltered <- sce.grun

This dataset lacks mitochondrial genes so we will do without them for quality control. We compute the median and MAD while blocking on the donor; for donors where the assumption of a majority of high-quality cells seems to be violated (Figure 5.1), we compute an appropriate threshold using the other donors as specified in the subset= argument.

library(scater)
stats <- perCellQCMetrics(sce.grun)

qc <- quickPerCellQC(stats, percent_subsets="altexps_ERCC_percent",
    batch=sce.grun$donor,
    subset=sce.grun$donor %in% c("D17", "D7", "D2"))

sce.grun <- sce.grun[,!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 of the Grun pancreas dataset. Each point represents a cell and is colored according to whether that cell was discarded.

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

colSums(as.matrix(qc), na.rm=TRUE)
##              low_lib_size            low_n_features high_altexps_ERCC_percent 
##                       450                       510                       606 
##                   discard 
##                       664

5.4 Normalization

library(scran)
set.seed(1000) # for irlba. 
clusters <- quickCluster(sce.grun)
sce.grun <- computeSumFactors(sce.grun, clusters=clusters)
sce.grun <- logNormCounts(sce.grun)
summary(sizeFactors(sce.grun))
##    Min. 1st Qu.  Median    Mean 3rd Qu.    Max. 
##  0.0933  0.5044  0.7894  1.0000  1.2307 10.8933
plot(librarySizeFactors(sce.grun), sizeFactors(sce.grun), pch=16,
    xlab="Library size factors", ylab="Deconvolution factors", log="xy")
Relationship between the library size factors and the deconvolution size factors in the Grun pancreas dataset.

Figure 5.2: Relationship between the library size factors and the deconvolution size factors in the Grun pancreas dataset.

5.5 Variance modelling

We block on a combined plate and donor factor.

block <- paste0(sce.grun$sample, "_", sce.grun$donor)
dec.grun <- modelGeneVarWithSpikes(sce.grun, spikes="ERCC", block=block)
top.grun <- getTopHVGs(dec.grun, prop=0.1)

We examine the number of cells in each level of the blocking factor.

table(block)
## block
##                  CD13+ sorted cells_D17       CD24+ CD44+ live sorted cells_D17 
##                                      87                                      87 
##                  CD63+ sorted cells_D10                TGFBR3+ sorted cells_D17 
##                                      40                                      90 
## exocrine fraction, live sorted cells_D2 exocrine fraction, live sorted cells_D3 
##                                      82                                       7 
##        live sorted cells, library 1_D10        live sorted cells, library 1_D17 
##                                      33                                      88 
##         live sorted cells, library 1_D3         live sorted cells, library 1_D7 
##                                      25                                      85 
##        live sorted cells, library 2_D10        live sorted cells, library 2_D17 
##                                      35                                      83 
##         live sorted cells, library 2_D3         live sorted cells, library 2_D7 
##                                      27                                      84 
##         live sorted cells, library 3_D3         live sorted cells, library 3_D7 
##                                      16                                      83 
##         live sorted cells, library 4_D3         live sorted cells, library 4_D7 
##                                      29                                      83
par(mfrow=c(6,3))
blocked.stats <- dec.grun$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 Grun 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 1.4: Per-gene variance as a function of the mean for the log-expression values in the Grun 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.

5.6 Data integration

library(batchelor)
set.seed(1001010)
merged.grun <- fastMNN(sce.grun, subset.row=top.grun, batch=sce.grun$donor)
metadata(merged.grun)$merge.info$lost.var
##           D10      D17       D2      D3      D7
## [1,] 0.030387 0.032160 0.000000 0.00000 0.00000
## [2,] 0.007869 0.013013 0.038314 0.00000 0.00000
## [3,] 0.004285 0.005418 0.008295 0.05409 0.00000
## [4,] 0.013813 0.016670 0.016313 0.01528 0.05643

5.7 Dimensionality reduction

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

5.8 Clustering

snn.gr <- buildSNNGraph(merged.grun, use.dimred="corrected")
colLabels(merged.grun) <- factor(igraph::cluster_walktrap(snn.gr)$membership)
table(Cluster=colLabels(merged.grun), Donor=merged.grun$batch)
##        Donor
## Cluster D10 D17  D2  D3  D7
##      1   32  71  31  80  29
##      2    3  10   3   3   7
##      3    1  10   0   0   8
##      4   11  70  28   2  69
##      5   11 119   0   0  55
##      6    3  42   0   0   9
##      7   16  37  15  11  45
##      8   14  30   3   2  66
##      9    5  18   0   2  34
##      10   5  13   0   0  10
##      11   3   2   2   4   2
##      12   4  13   0   0   1
gridExtra::grid.arrange(
    plotTSNE(merged.grun, colour_by="label"),
    plotTSNE(merged.grun, colour_by="batch"),
    ncol=2
)
Obligatory $t$-SNE plots of the Grun pancreas dataset. Each point represents a cell that is colored by cluster (left) or batch (right).

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

Session Info

R version 4.6.0 RC (2026-04-17 r89917)
Platform: x86_64-pc-linux-gnu
Running under: Ubuntu 24.04.4 LTS

Matrix products: default
BLAS:   /home/biocbuild/bbs-3.24-bioc/R/lib/libRblas.so 
LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.12.0  LAPACK version 3.12.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] batchelor_1.29.0            scran_1.41.0               
 [3] scater_1.41.1               ggplot2_4.0.3              
 [5] scuttle_1.23.0              org.Hs.eg.db_3.23.1        
 [7] AnnotationDbi_1.75.0        scRNAseq_2.27.0            
 [9] SingleCellExperiment_1.35.0 SummarizedExperiment_1.43.0
[11] Biobase_2.73.1              GenomicRanges_1.65.0       
[13] Seqinfo_1.3.0               IRanges_2.47.0             
[15] S4Vectors_0.51.1            BiocGenerics_0.59.0        
[17] generics_0.1.4              MatrixGenerics_1.25.0      
[19] matrixStats_1.5.0           BiocStyle_2.41.0           
[21] rebook_1.23.0              

loaded via a namespace (and not attached):
  [1] RColorBrewer_1.1-3        jsonlite_2.0.0           
  [3] CodeDepends_0.6.7         magrittr_2.0.5           
  [5] ggbeeswarm_0.7.3          GenomicFeatures_1.65.0   
  [7] gypsum_1.9.0              farver_2.1.2             
  [9] rmarkdown_2.31            BiocIO_1.23.3            
 [11] vctrs_0.7.3               DelayedMatrixStats_1.35.0
 [13] memoise_2.0.1             Rsamtools_2.29.0         
 [15] RCurl_1.98-1.18           htmltools_0.5.9          
 [17] S4Arrays_1.13.0           BiocBaseUtils_1.15.0     
 [19] AnnotationHub_4.3.0       curl_7.1.0               
 [21] BiocNeighbors_2.7.0       Rhdf5lib_2.1.0           
 [23] SparseArray_1.13.2        rhdf5_2.57.0             
 [25] sass_0.4.10               alabaster.base_1.13.0    
 [27] bslib_0.10.0              alabaster.sce_1.13.0     
 [29] httr2_1.2.2               cachem_1.1.0             
 [31] ResidualMatrix_1.23.0     GenomicAlignments_1.49.0 
 [33] igraph_2.3.1              lifecycle_1.0.5          
 [35] pkgconfig_2.0.3           rsvd_1.0.5               
 [37] Matrix_1.7-5              R6_2.6.1                 
 [39] fastmap_1.2.0             digest_0.6.39            
 [41] dqrng_0.4.1               irlba_2.3.7              
 [43] ExperimentHub_3.3.0       RSQLite_2.4.6            
 [45] beachmat_2.29.0           labeling_0.4.3           
 [47] filelock_1.0.3            httr_1.4.8               
 [49] abind_1.4-8               compiler_4.6.0           
 [51] bit64_4.8.0               withr_3.0.2              
 [53] S7_0.2.2                  BiocParallel_1.47.0      
 [55] viridis_0.6.5             DBI_1.3.0                
 [57] HDF5Array_1.41.0          alabaster.ranges_1.13.0  
 [59] alabaster.schemas_1.13.0  rappdirs_0.3.4           
 [61] DelayedArray_0.39.1       bluster_1.23.0           
 [63] rjson_0.2.23              tools_4.6.0              
 [65] vipor_0.4.7               otel_0.2.0               
 [67] beeswarm_0.4.0            glue_1.8.1               
 [69] h5mread_1.5.0             restfulr_0.0.16          
 [71] rhdf5filters_1.25.0       grid_4.6.0               
 [73] Rtsne_0.17                cluster_2.1.8.2          
 [75] gtable_0.3.6              ensembldb_2.37.0         
 [77] metapod_1.21.0            BiocSingular_1.29.0      
 [79] ScaledMatrix_1.21.0       XVector_0.53.0           
 [81] ggrepel_0.9.8             BiocVersion_3.24.0       
 [83] pillar_1.11.1             limma_3.69.0             
 [85] dplyr_1.2.1               BiocFileCache_3.3.0      
 [87] lattice_0.22-9            rtracklayer_1.73.0       
 [89] bit_4.6.0                 tidyselect_1.2.1         
 [91] locfit_1.5-9.12           Biostrings_2.81.1        
 [93] knitr_1.51                gridExtra_2.3            
 [95] bookdown_0.46             ProtGenerics_1.45.0      
 [97] edgeR_4.11.0              xfun_0.57                
 [99] statmod_1.5.1             UCSC.utils_1.9.0         
[101] lazyeval_0.2.3            yaml_2.3.12              
[103] evaluate_1.0.5            codetools_0.2-20         
[105] cigarillo_1.3.0           tibble_3.3.1             
[107] alabaster.matrix_1.13.0   BiocManager_1.30.27      
[109] graph_1.91.0              cli_3.6.6                
[111] jquerylib_0.1.4           dichromat_2.0-0.1        
[113] Rcpp_1.1.1-1.1            GenomeInfoDb_1.49.0      
[115] dir.expiry_1.21.0         dbplyr_2.5.2             
[117] png_0.1-9                 XML_3.99-0.23            
[119] parallel_4.6.0            blob_1.3.0               
[121] AnnotationFilter_1.37.0   sparseMatrixStats_1.25.0 
[123] bitops_1.0-9              viridisLite_0.4.3        
[125] alabaster.se_1.13.0       scales_1.4.0             
[127] crayon_1.5.3              rlang_1.2.0              
[129] cowplot_1.2.0             KEGGREST_1.53.0          

References

Grun, D., M. J. Muraro, J. C. Boisset, K. Wiebrands, A. Lyubimova, G. Dharmadhikari, M. van den Born, et al. 2016. β€œDe Novo Prediction of Stem Cell Identity using Single-Cell Transcriptome Data.” Cell Stem Cell 19 (2): 266–77.