In this vignette we will demonstrate how to construct cell trajectories with Monocle 3 using single-cell ATAC-seq data. Please see the Monocle 3 website for information about installing Monocle 3.

To facilitate conversion between the Seurat (used by Signac) and CellDataSet (used by Monocle 3) formats, we will use a conversion function in the SeuratWrappers package available on GitHub.

Data loading

We will use a single-cell ATAC-seq dataset containing human CD34+ hematopoietic stem and progenitor cells published by Satpathy and Granja et al. (2019, Nature Biotechnology).

The processed dataset is available on NCBI GEO here: https://www.ncbi.nlm.nih.gov/geo/query/acc.cgi?acc=GSE129785

Note that the fragment file is present inside the GSE129785_RAW.tar archive, and the index for the fragment file is not supplied. You can index the file yourself using tabix, for example: tabix -p bed <fragment_file>.

First we will load the dataset and perform some standard preprocessing using Signac.

filepath <- "~/data/satpathy19/GSE129785/GSE129785_scATAC-Hematopoiesis-CD34"

peaks <- read.table(paste0(filepath, ".peaks.txt.gz"), header = TRUE)
cells <- read.table(paste0(filepath, ".cell_barcodes.txt.gz"), header = TRUE, stringsAsFactors = FALSE)
rownames(cells) <- make.unique(cells$Barcodes)

mtx <- readMM(file = paste0(filepath, ".mtx"))
mtx <- as(object = mtx, Class = "dgCMatrix")
colnames(mtx) <- rownames(cells)
rownames(mtx) <- peaks$Feature
bone_assay <- CreateChromatinAssay(
  counts = mtx,
  min.cells = 5,
  fragments = "~/data/satpathy19/GSE129785/GSM3722029_CD34_Progenitors_Rep1_fragments.tsv.gz",
  sep = c("_", "_"),
  genome = "hg19"
)
bone <- CreateSeuratObject(
  counts = bone_assay,
  meta.data = cells,
  assay = "ATAC"
)

# The dataset contains multiple cell types
# We can subset to include just one replicate of CD34+ progenitor cells
bone <- bone[, bone$Group_Barcode == "CD34_Progenitors_Rep1"]

# add cell type annotations from the original paper
cluster_names <- c("HSC",   "MEP",  "CMP-BMP",  "LMPP", "CLP",  "Pro-B",    "Pre-B",    "GMP",
                  "MDP",    "pDC",  "cDC",  "Monocyte-1",   "Monocyte-2",   "Naive-B",  "Memory-B",
                  "Plasma-cell",    "Basophil", "Immature-NK",  "Mature-NK1",   "Mature-NK2",   "Naive-CD4-T1",
                  "Naive-CD4-T2",   "Naive-Treg",   "Memory-CD4-T", "Treg", "Naive-CD8-T1", "Naive-CD8-T2",
                  "Naive-CD8-T3",   "Central-memory-CD8-T", "Effector-memory-CD8-T",    "Gamma delta T")
num.labels <- length(cluster_names)
names(cluster_names) <- paste0( rep("Cluster", num.labels), seq(num.labels) )
bone$celltype <- cluster_names[as.character(bone$Clusters)]

bone[["ATAC"]]
## ChromatinAssay data with 570650 features for 9913 cells
## Variable features: 0 
## Genome: hg19 
## Annotation present: FALSE 
## Motifs present: FALSE 
## Fragment files: 1

Next we can add gene annotations for the hg19 genome to the object. This will be useful for computing quality control metrics (TSS enrichment score) and plotting.

library(EnsDb.Hsapiens.v75)

# extract gene annotations from EnsDb
annotations <- GetGRangesFromEnsDb(ensdb = EnsDb.Hsapiens.v75)

# change to UCSC style since the data was mapped to hg19
seqlevelsStyle(annotations) <- 'UCSC'
genome(annotations) <- "hg19"

# add the gene information to the object
Annotation(bone) <- annotations

Quality control

We’ll compute TSS enrichment, nucleosome signal score, and the percentage of counts in genomic blacklist regions for each cell, and use these metrics to help remove low quality cells from the datasets.

bone <- TSSEnrichment(bone)
bone <- NucleosomeSignal(bone)
bone$blacklist_fraction <- FractionCountsInRegion(bone, regions = blacklist_hg19)
VlnPlot(
  object = bone,
  features = c("nCount_ATAC", "TSS.enrichment", "nucleosome_signal", "blacklist_fraction"),
  pt.size = 0.1,
  ncol = 4
)

bone <- bone[, (bone$nCount_ATAC < 50000) &
               (bone$TSS.enrichment > 2) & 
               (bone$nucleosome_signal < 5)]

Dataset preprocessing

Next we can run a standard scATAC-seq analysis pipeline using Signac to perform dimension reduction, clustering, and cell type annotation.

bone <- FindTopFeatures(bone, min.cells = 10)
bone <- RunTFIDF(bone)
bone <- RunSVD(bone, n = 100)
DepthCor(bone)

bone <- RunUMAP(
  bone,
  reduction = "lsi",
  dims = 2:50,
  reduction.name = "UMAP"
)
bone <- FindNeighbors(bone, dims = 2:50, reduction = "lsi")
bone <- FindClusters(bone, resolution = 0.8, algorithm = 3)
## Modularity Optimizer version 1.3.0 by Ludo Waltman and Nees Jan van Eck
## 
## Number of nodes: 9765
## Number of edges: 433498
## 
## Running smart local moving algorithm...
## Maximum modularity in 10 random starts: 0.8474
## Number of communities: 18
## Elapsed time: 5 seconds
DimPlot(bone, label = TRUE) + NoLegend()

Assign each cluster to the most common cell type based on the original annotations from the paper.

for(i in levels(bone)) {
  cells_to_reid <- WhichCells(bone, idents = i)
  newid <- names(sort(table(bone$celltype[cells_to_reid]),decreasing=TRUE))[1]
  Idents(bone, cells = cells_to_reid) <- newid
}
bone$assigned_celltype <- Idents(bone)
DimPlot(bone, label = TRUE)

Next we can subset the different lineages and create a trajectory for each lineage. Another way to build the trajectories is to use the whole dataset and build separate pseudotime trajectories for the different cell partitions found by Monocle 3.

DefaultAssay(bone) <- "ATAC"

erythroid <- bone[,  bone$assigned_celltype %in% c("HSC", "MEP", "CMP-BMP")]
lymphoid <- bone[, bone$assigned_celltype %in% c("HSC", "LMPP", "GMP", "CLP", "Pro-B", "pDC", "MDP", "GMP")]

Building trajectories with Monocle 3

We can convert the Seurat object to a CellDataSet object using the as.cell_data_set() function from SeuratWrappers and build the trajectories using Monocle 3. We’ll do this separately for erythroid and lymphoid lineages, but you could explore other strategies building a trajectory for all lineages together.

erythroid.cds <- as.cell_data_set(erythroid)
erythroid.cds <- cluster_cells(cds = erythroid.cds, reduction_method = "UMAP")
erythroid.cds <- learn_graph(erythroid.cds, use_partition = TRUE)

lymphoid.cds <- as.cell_data_set(lymphoid)
lymphoid.cds <- cluster_cells(cds = lymphoid.cds, reduction_method = "UMAP")
lymphoid.cds <- learn_graph(lymphoid.cds, use_partition = TRUE)

To compute pseudotime estimates for each trajectory we need to decide what the start of each trajectory is. In our case, we know that the hematopoietic stem cells are the progenitors of other cell types in the trajectory, so we can set these cells as the root of the trajectory. Monocle 3 includes an interactive function to select cells as the root nodes in the graph. This function will be launched if calling order_cells() without specifying the root_cells parameter. Here we’ve pre-selected some cells as the root, and saved these to a file for reproducibility. This file can be downloaded here.

# load the pre-selected HSCs
hsc <- readLines("../vignette_data/hsc_cells.txt")
# order cells
erythroid.cds <- order_cells(erythroid.cds, reduction_method = "UMAP", root_cells = hsc)
lymphoid.cds <- order_cells(lymphoid.cds, reduction_method = "UMAP", root_cells = hsc)

# plot trajectories colored by pseudotime
plot_cells(
  cds = erythroid.cds,
  color_cells_by = "pseudotime",
  show_trajectory_graph = TRUE
)

plot_cells(
  cds = lymphoid.cds,
  color_cells_by = "pseudotime",
  show_trajectory_graph = TRUE
)

Extract the pseudotime values and add to the Seurat object

bone <- AddMetaData(
  object = bone,
  metadata = erythroid.cds@principal_graph_aux@listData$UMAP$pseudotime,
  col.name = "Erythroid"
)

bone <- AddMetaData(
  object = bone,
  metadata = lymphoid.cds@principal_graph_aux@listData$UMAP$pseudotime,
  col.name = "Lymphoid"
)
FeaturePlot(bone, c("Erythroid", "Lymphoid"), pt.size = 0.1) & scale_color_viridis_c()

Acknowledgements

Thanks to the developers of Monocle 3, especially Cole Trapnell, Hannah Pliner, and members of the Trapnell lab. If you use Monocle please cite the Monocle papers.

Session Info
## R version 4.0.1 (2020-06-06)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 18.04.5 LTS
## 
## Matrix products: default
## BLAS:   /usr/lib/x86_64-linux-gnu/blas/libblas.so.3.7.1
## LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.7.1
## 
## 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       
## 
## attached base packages:
## [1] stats4    parallel  stats     graphics  grDevices utils     datasets 
## [8] methods   base     
## 
## other attached packages:
##  [1] EnsDb.Hsapiens.v75_2.99.0   ensembldb_2.14.0           
##  [3] AnnotationFilter_1.14.0     GenomicFeatures_1.42.3     
##  [5] AnnotationDbi_1.52.0        patchwork_1.1.1            
##  [7] ggplot2_3.3.3               Matrix_1.3-2               
##  [9] monocle3_1.0.0              SingleCellExperiment_1.12.0
## [11] SummarizedExperiment_1.20.0 GenomicRanges_1.42.0       
## [13] GenomeInfoDb_1.26.5         IRanges_2.24.1             
## [15] S4Vectors_0.28.1            MatrixGenerics_1.2.1       
## [17] matrixStats_0.58.0          Biobase_2.50.0             
## [19] BiocGenerics_0.36.0         SeuratWrappers_0.3.0       
## [21] SeuratObject_4.0.0          Seurat_4.0.1.9005          
## [23] Signac_1.2.0               
## 
## loaded via a namespace (and not attached):
##   [1] rappdirs_0.3.3           SnowballC_0.7.0          rtracklayer_1.50.0      
##   [4] scattermore_0.7          ragg_1.1.2               tidyr_1.1.3             
##   [7] bit64_4.0.5              knitr_1.33               irlba_2.3.3             
##  [10] DelayedArray_0.16.3      data.table_1.14.0        rpart_4.1-15            
##  [13] RCurl_1.98-1.3           generics_0.1.0           leidenbase_0.1.2        
##  [16] cowplot_1.1.1            RSQLite_2.2.7            RANN_2.6.1              
##  [19] proxy_0.4-25             future_1.21.0            bit_4.0.4               
##  [22] spatstat.data_2.1-0      xml2_1.3.2               httpuv_1.6.0            
##  [25] assertthat_0.2.1         viridis_0.6.0            xfun_0.22               
##  [28] hms_1.0.0                jquerylib_0.1.4          evaluate_0.14           
##  [31] promises_1.2.0.1         fansi_0.4.2              progress_1.2.2          
##  [34] dbplyr_2.1.1             igraph_1.2.6             DBI_1.1.1               
##  [37] htmlwidgets_1.5.3        sparsesvd_0.2            spatstat.geom_2.1-0     
##  [40] purrr_0.3.4              ellipsis_0.3.1           RSpectra_0.16-0         
##  [43] dplyr_1.0.5              backports_1.2.1          biomaRt_2.46.3          
##  [46] deldir_0.2-10            vctrs_0.3.7              remotes_2.3.0           
##  [49] ROCR_1.0-11              abind_1.4-5              cachem_1.0.4            
##  [52] withr_2.4.2              ggforce_0.3.3            BSgenome_1.58.0         
##  [55] checkmate_2.0.0          sctransform_0.3.2        GenomicAlignments_1.26.0
##  [58] prettyunits_1.1.1        goftest_1.2-2            cluster_2.1.2           
##  [61] lazyeval_0.2.2           crayon_1.4.1             pkgconfig_2.0.3         
##  [64] slam_0.1-48              labeling_0.4.2           tweenr_1.0.2            
##  [67] nlme_3.1-152             ProtGenerics_1.22.0      nnet_7.3-15             
##  [70] rlang_0.4.10             globals_0.14.0           lifecycle_1.0.0         
##  [73] miniUI_0.1.1.1           BiocFileCache_1.14.0     rsvd_1.0.5              
##  [76] dichromat_2.0-0          rprojroot_2.0.2          polyclip_1.10-0         
##  [79] lmtest_0.9-38            ggseqlogo_0.1            zoo_1.8-9               
##  [82] base64enc_0.1-3          ggridges_0.5.3           png_0.1-7               
##  [85] viridisLite_0.4.0        bitops_1.0-7             KernSmooth_2.23-18      
##  [88] Biostrings_2.58.0        blob_1.2.1               stringr_1.4.0           
##  [91] parallelly_1.24.0        jpeg_0.1-8.1             scales_1.1.1            
##  [94] memoise_2.0.0            magrittr_2.0.1           plyr_1.8.6              
##  [97] ica_1.0-2                zlibbioc_1.36.0          compiler_4.0.1          
## [100] RColorBrewer_1.1-2       fitdistrplus_1.1-3       Rsamtools_2.6.0         
## [103] XVector_0.30.0           listenv_0.8.0            pbapply_1.4-3           
## [106] htmlTable_2.1.0          Formula_1.2-4            MASS_7.3-53.1           
## [109] mgcv_1.8-33              tidyselect_1.1.0         stringi_1.5.3           
## [112] textshaping_0.3.3        highr_0.9                yaml_2.2.1              
## [115] askpass_1.1              latticeExtra_0.6-29      ggrepel_0.9.1           
## [118] grid_4.0.1               sass_0.3.1               VariantAnnotation_1.36.0
## [121] fastmatch_1.1-0          tools_4.0.1              future.apply_1.7.0      
## [124] rstudioapi_0.13          foreign_0.8-81           lsa_0.73.2              
## [127] gridExtra_2.3            farver_2.1.0             Rtsne_0.15              
## [130] digest_0.6.27            BiocManager_1.30.12      shiny_1.6.0             
## [133] qlcMatrix_0.9.7          Rcpp_1.0.6               later_1.2.0             
## [136] RcppAnnoy_0.0.18         httr_1.4.2               biovizBase_1.38.0       
## [139] colorspace_2.0-0         XML_3.99-0.6             fs_1.5.0                
## [142] tensor_1.5               reticulate_1.19          splines_4.0.1           
## [145] uwot_0.1.10              RcppRoll_0.3.0           spatstat.utils_2.1-0    
## [148] pkgdown_1.6.1            plotly_4.9.3             systemfonts_1.0.1       
## [151] xtable_1.8-4             jsonlite_1.7.2           R6_2.5.0                
## [154] Hmisc_4.5-0              pillar_1.6.0             htmltools_0.5.1.1       
## [157] mime_0.10                glue_1.4.2               fastmap_1.1.0           
## [160] BiocParallel_1.24.1      codetools_0.2-18         utf8_1.2.1              
## [163] lattice_0.20-41          bslib_0.2.4              spatstat.sparse_2.0-0   
## [166] tibble_3.1.1             curl_4.3                 leiden_0.3.7            
## [169] openssl_1.4.3            survival_3.2-11          rmarkdown_2.7           
## [172] docopt_0.7.1             desc_1.3.0               munsell_0.5.0           
## [175] GenomeInfoDbData_1.2.4   reshape2_1.4.4           gtable_0.3.0            
## [178] spatstat.core_2.1-2