Chapter 11 Paul mouse HSC (MARS-seq)

11.1 Introduction

This performs an analysis of the mouse haematopoietic stem cell (HSC) dataset generated with MARS-seq (Paul et al. 2015). Cells were extracted from multiple mice under different experimental conditions (i.e., sorting protocols) and libraries were prepared using a series of 384-well plates.

11.2 Data loading

library(scRNAseq)
sce.paul <- PaulHSCData(ensembl=TRUE)
library(AnnotationHub)
ens.mm.v97 <- AnnotationHub()[["AH73905"]]
anno <- select(ens.mm.v97, keys=rownames(sce.paul), 
    keytype="GENEID", columns=c("SYMBOL", "SEQNAME"))
rowData(sce.paul) <- anno[match(rownames(sce.paul), anno$GENEID),]

After loading and annotation, we inspect the resulting SingleCellExperiment object:

sce.paul
## class: SingleCellExperiment 
## dim: 17483 10368 
## metadata(0):
## assays(1): counts
## rownames(17483): ENSMUSG00000007777 ENSMUSG00000107002 ...
##   ENSMUSG00000039068 ENSMUSG00000064363
## rowData names(3): GENEID SYMBOL SEQNAME
## colnames(10368): W29953 W29954 ... W76335 W76336
## colData names(13): Well_ID Seq_batch_ID ... CD34_measurement
##   FcgR3_measurement
## reducedDimNames(0):
## mainExpName: NULL
## altExpNames(0):

11.3 Quality control

unfiltered <- sce.paul

For some reason, only one mitochondrial transcripts are available, so we will perform quality control using only the library size and number of detected features. Ideally, we would simply block on the plate of origin to account for differences in processing, but unfortunately, it seems that many plates have a large proportion (if not outright majority) of cells with poor values for both metrics. We identify such plates based on the presence of very low outlier thresholds, for some arbitrary definition of “low”; we then redefine thresholds using information from the other (presumably high-quality) plates.

library(scater)
stats <- perCellQCMetrics(sce.paul) 
qc <- quickPerCellQC(stats, batch=sce.paul$Plate_ID)

# Detecting batches with unusually low threshold values.
lib.thresholds <- attr(qc$low_lib_size, "thresholds")["lower",]
nfeat.thresholds <- attr(qc$low_n_features, "thresholds")["lower",]
ignore <- union(names(lib.thresholds)[lib.thresholds < 100],
    names(nfeat.thresholds)[nfeat.thresholds < 100])

# Repeating the QC using only the "high-quality" batches.
qc2 <- quickPerCellQC(stats, batch=sce.paul$Plate_ID,
    subset=!sce.paul$Plate_ID %in% ignore)
sce.paul <- sce.paul[,!qc2$discard]

We examine the number of cells discarded for each reason.

colSums(as.matrix(qc2))
##   low_lib_size low_n_features        discard 
##           1695           1781           1783

We create some diagnostic plots for each metric (Figure 11.1).

colData(unfiltered) <- cbind(colData(unfiltered), stats)
unfiltered$discard <- qc2$discard
unfiltered$Plate_ID <- factor(unfiltered$Plate_ID)

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

Figure 11.1: Distribution of each QC metric across cells in the Paul HSC dataset. Each point represents a cell and is colored according to whether that cell was discarded.

11.4 Normalization

library(scran)
set.seed(101000110)
clusters <- quickCluster(sce.paul)
sce.paul <- computeSumFactors(sce.paul, clusters=clusters)
sce.paul <- logNormCounts(sce.paul)

We examine some key metrics for the distribution of size factors, and compare it to the library sizes as a sanity check (Figure 11.2).

summary(sizeFactors(sce.paul))
##    Min. 1st Qu.  Median    Mean 3rd Qu.    Max. 
##   0.057   0.422   0.775   1.000   1.335   9.654
plot(librarySizeFactors(sce.paul), sizeFactors(sce.paul), pch=16,
    xlab="Library size factors", ylab="Deconvolution factors", log="xy")
Relationship between the library size factors and the deconvolution size factors in the Paul HSC dataset.

Figure 11.2: Relationship between the library size factors and the deconvolution size factors in the Paul HSC dataset.

11.5 Variance modelling

We fit a mean-variance trend to the endogenous genes to detect highly variable genes. Unfortunately, the plates are confounded with an experimental treatment (Batch_desc) so we cannot block on the plate of origin.

set.seed(00010101)
dec.paul <- modelGeneVarByPoisson(sce.paul)
top.paul <- getTopHVGs(dec.paul, prop=0.1)
plot(dec.paul$mean, dec.paul$total, pch=16, cex=0.5, 
    xlab="Mean of log-expression", ylab="Variance of log-expression")
curve(metadata(dec.paul)$trend(x), col="blue", add=TRUE)
Per-gene variance as a function of the mean for the log-expression values in the Paul HSC dataset. Each point represents a gene (black) with the mean-variance trend (blue) fitted to simulated Poisson noise.

Figure 11.3: Per-gene variance as a function of the mean for the log-expression values in the Paul HSC dataset. Each point represents a gene (black) with the mean-variance trend (blue) fitted to simulated Poisson noise.

11.6 Dimensionality reduction

set.seed(101010011)
sce.paul <- denoisePCA(sce.paul, technical=dec.paul, subset.row=top.paul)
sce.paul <- runTSNE(sce.paul, dimred="PCA")

We check that the number of retained PCs is sensible.

ncol(reducedDim(sce.paul, "PCA"))
## [1] 13

11.7 Clustering

snn.gr <- buildSNNGraph(sce.paul, use.dimred="PCA", type="jaccard")
colLabels(sce.paul) <- factor(igraph::cluster_louvain(snn.gr)$membership)

These is a strong relationship between the cluster and the experimental treatment (Figure 11.4), which is to be expected. Of course, this may also be attributable to some batch effect; the confounded nature of the experimental design makes it difficult to make any confident statements either way.

tab <- table(colLabels(sce.paul), sce.paul$Batch_desc)
rownames(tab) <- paste("Cluster", rownames(tab))
pheatmap::pheatmap(log10(tab+10), color=viridis::viridis(100))
Heatmap of the distribution of cells across clusters (rows) for each experimental treatment (column).

Figure 11.4: Heatmap of the distribution of cells across clusters (rows) for each experimental treatment (column).

plotTSNE(sce.paul, colour_by="label")
Obligatory $t$-SNE plot of the Paul HSC dataset, where each point represents a cell and is colored according to the assigned cluster.

Figure 11.5: Obligatory \(t\)-SNE plot of the Paul HSC dataset, where each point represents a cell and is colored according to the assigned cluster.

plotTSNE(sce.paul, colour_by="label", other_fields="Batch_desc") + 
    facet_wrap(~Batch_desc)
Obligatory $t$-SNE plot of the Paul HSC dataset faceted by the treatment condition, where each point represents a cell and is colored according to the assigned cluster.

Figure 11.6: Obligatory \(t\)-SNE plot of the Paul HSC dataset faceted by the treatment condition, where each point represents a cell and is colored according to the assigned cluster.

Session Info

R version 4.3.0 RC (2023-04-13 r84269)
Platform: x86_64-pc-linux-gnu (64-bit)
Running under: Ubuntu 22.04.2 LTS

Matrix products: default
BLAS:   /home/biocbuild/bbs-3.17-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] scran_1.28.0                scater_1.28.0              
 [3] ggplot2_3.4.2               scuttle_1.10.0             
 [5] AnnotationHub_3.8.0         BiocFileCache_2.8.0        
 [7] dbplyr_2.3.2                ensembldb_2.24.0           
 [9] AnnotationFilter_1.24.0     GenomicFeatures_1.52.0     
[11] AnnotationDbi_1.62.0        scRNAseq_2.13.0            
[13] SingleCellExperiment_1.22.0 SummarizedExperiment_1.30.0
[15] Biobase_2.60.0              GenomicRanges_1.52.0       
[17] GenomeInfoDb_1.36.0         IRanges_2.34.0             
[19] S4Vectors_0.38.0            BiocGenerics_0.46.0        
[21] MatrixGenerics_1.12.0       matrixStats_0.63.0         
[23] BiocStyle_2.28.0            rebook_1.10.0              

loaded via a namespace (and not attached):
  [1] RColorBrewer_1.1-3            jsonlite_1.8.4               
  [3] CodeDepends_0.6.5             magrittr_2.0.3               
  [5] ggbeeswarm_0.7.1              farver_2.1.1                 
  [7] rmarkdown_2.21                BiocIO_1.10.0                
  [9] zlibbioc_1.46.0               vctrs_0.6.2                  
 [11] memoise_2.0.1                 Rsamtools_2.16.0             
 [13] DelayedMatrixStats_1.22.0     RCurl_1.98-1.12              
 [15] htmltools_0.5.5               progress_1.2.2               
 [17] curl_5.0.0                    BiocNeighbors_1.18.0         
 [19] sass_0.4.5                    bslib_0.4.2                  
 [21] cachem_1.0.7                  GenomicAlignments_1.36.0     
 [23] igraph_1.4.2                  mime_0.12                    
 [25] lifecycle_1.0.3               pkgconfig_2.0.3              
 [27] rsvd_1.0.5                    Matrix_1.5-4                 
 [29] R6_2.5.1                      fastmap_1.1.1                
 [31] GenomeInfoDbData_1.2.10       shiny_1.7.4                  
 [33] digest_0.6.31                 colorspace_2.1-0             
 [35] dqrng_0.3.0                   irlba_2.3.5.1                
 [37] ExperimentHub_2.8.0           RSQLite_2.3.1                
 [39] beachmat_2.16.0               labeling_0.4.2               
 [41] filelock_1.0.2                fansi_1.0.4                  
 [43] httr_1.4.5                    compiler_4.3.0               
 [45] bit64_4.0.5                   withr_2.5.0                  
 [47] BiocParallel_1.34.0           viridis_0.6.2                
 [49] DBI_1.1.3                     highr_0.10                   
 [51] biomaRt_2.56.0                rappdirs_0.3.3               
 [53] DelayedArray_0.26.0           bluster_1.10.0               
 [55] rjson_0.2.21                  tools_4.3.0                  
 [57] vipor_0.4.5                   beeswarm_0.4.0               
 [59] interactiveDisplayBase_1.38.0 httpuv_1.6.9                 
 [61] glue_1.6.2                    restfulr_0.0.15              
 [63] promises_1.2.0.1              grid_4.3.0                   
 [65] Rtsne_0.16                    cluster_2.1.4                
 [67] generics_0.1.3                gtable_0.3.3                 
 [69] hms_1.1.3                     metapod_1.8.0                
 [71] BiocSingular_1.16.0           ScaledMatrix_1.8.0           
 [73] xml2_1.3.3                    utf8_1.2.3                   
 [75] XVector_0.40.0                ggrepel_0.9.3                
 [77] BiocVersion_3.17.1            pillar_1.9.0                 
 [79] stringr_1.5.0                 limma_3.56.0                 
 [81] later_1.3.0                   dplyr_1.1.2                  
 [83] lattice_0.21-8                rtracklayer_1.60.0           
 [85] bit_4.0.5                     tidyselect_1.2.0             
 [87] locfit_1.5-9.7                Biostrings_2.68.0            
 [89] knitr_1.42                    gridExtra_2.3                
 [91] bookdown_0.33                 ProtGenerics_1.32.0          
 [93] edgeR_3.42.0                  xfun_0.39                    
 [95] statmod_1.5.0                 pheatmap_1.0.12              
 [97] stringi_1.7.12                lazyeval_0.2.2               
 [99] yaml_2.3.7                    evaluate_0.20                
[101] codetools_0.2-19              tibble_3.2.1                 
[103] BiocManager_1.30.20           graph_1.78.0                 
[105] cli_3.6.1                     xtable_1.8-4                 
[107] munsell_0.5.0                 jquerylib_0.1.4              
[109] Rcpp_1.0.10                   dir.expiry_1.8.0             
[111] png_0.1-8                     XML_3.99-0.14                
[113] parallel_4.3.0                ellipsis_0.3.2               
[115] blob_1.2.4                    prettyunits_1.1.1            
[117] sparseMatrixStats_1.12.0      bitops_1.0-7                 
[119] viridisLite_0.4.1             scales_1.2.1                 
[121] purrr_1.0.1                   crayon_1.5.2                 
[123] rlang_1.1.0                   cowplot_1.1.1                
[125] KEGGREST_1.40.0              

References

Paul, F., Y. Arkin, A. Giladi, D. A. Jaitin, E. Kenigsberg, H. Keren-Shaul, D. Winter, et al. 2015. “Transcriptional Heterogeneity and Lineage Commitment in Myeloid Progenitors.” Cell 163 (7): 1663–77.