vignettes/mouse_brain_vignette.Rmd
mouse_brain_vignette.Rmd
For this tutorial, we will be analyzing a single-cell ATAC-seq dataset of adult mouse brain cells provided by 10x Genomics. The following files are used in this vignette, all available through the 10x Genomics website.
View data download code
To download the required data, run the following lines in a shell:
wget http://cf.10xgenomics.com/samples/cell-atac/1.1.0/atac_v1_adult_brain_fresh_5k/atac_v1_adult_brain_fresh_5k_filtered_peak_bc_matrix.h5
wget http://cf.10xgenomics.com/samples/cell-atac/1.1.0/atac_v1_adult_brain_fresh_5k/atac_v1_adult_brain_fresh_5k_singlecell.csv
wget http://cf.10xgenomics.com/samples/cell-atac/1.1.0/atac_v1_adult_brain_fresh_5k/atac_v1_adult_brain_fresh_5k_fragments.tsv.gz
wget http://cf.10xgenomics.com/samples/cell-atac/1.1.0/atac_v1_adult_brain_fresh_5k/atac_v1_adult_brain_fresh_5k_fragments.tsv.gz.tbi
This vignette echoes the commands run in the introductory Signac vignette on human PBMC. We provide the same analysis in a different system to demonstrate performance and applicability to other tissue types, and to provide an example from another species.
First load in Signac, Seurat, and some other packages we will be using for analyzing mouse data.
library(Signac)
library(Seurat)
library(GenomeInfoDb)
library(EnsDb.Mmusculus.v79)
library(ggplot2)
library(patchwork)
set.seed(1234)
counts <- Read10X_h5("../vignette_data/atac_v1_adult_brain_fresh_5k_filtered_peak_bc_matrix.h5")
## Warning in sparseMatrix(i = indices[] + 1, p = indptr[], x = as.numeric(x =
## counts[]), : 'giveCsparse' has been deprecated; setting 'repr = "T"' for you
metadata <- read.csv(
file = "../vignette_data/atac_v1_adult_brain_fresh_5k_singlecell.csv",
header = TRUE,
row.names = 1
)
brain_assay <- CreateChromatinAssay(
counts = counts,
sep = c(":", "-"),
genome = "mm10",
fragments = '../vignette_data/atac_v1_adult_brain_fresh_5k_fragments.tsv.gz',
min.cells = 1
)
## Computing hash
brain <- CreateSeuratObject(
counts = brain_assay,
assay = 'peaks',
project = 'ATAC',
meta.data = metadata
)
## Warning in CreateSeuratObject.Assay(counts = brain_assay, assay = "peaks", :
## Some cells in meta.data not present in provided counts matrix.
## Warning: Keys should be one or more alphanumeric characters followed by an
## underscore, setting key from peaks to peaks_
We can also add gene annotations to the brain
object for the mouse genome. This will allow downstream functions to pull the gene annotation information directly from the object.
# extract gene annotations from EnsDb
annotations <- GetGRangesFromEnsDb(ensdb = EnsDb.Mmusculus.v79)
# change to UCSC style since the data was mapped to hg19
seqlevelsStyle(annotations) <- 'UCSC'
# add the gene information to the object
Annotation(brain) <- annotations
Next we compute some useful per-cell QC metrics.
brain <- NucleosomeSignal(object = brain)
We can look at the fragment length periodicity for all the cells, and group by cells with high or low nucleosomal signal strength. You can see that cells which are outliers for the mononucleosomal/ nucleosome-free ratio have different banding patterns. The remaining cells exhibit a pattern that is typical for a successful ATAC-seq experiment.
brain$nucleosome_group <- ifelse(brain$nucleosome_signal > 4, 'NS > 4', 'NS < 4')
FragmentHistogram(object = brain, group.by = 'nucleosome_group', region = 'chr1-1-10000000')
The enrichment of Tn5 integration events at transcriptional start sites (TSSs) can also be an important quality control metric to assess the targeting of Tn5 in ATAC-seq experiments. The ENCODE consortium defined a TSS enrichment score as the number of Tn5 integration site around the TSS normalized to the number of Tn5 integration sites in flanking regions. See the ENCODE documentation for more information about the TSS enrichment score (https://www.encodeproject.org/data-standards/terms/). We can calculate the TSS enrichment score for each cell using the TSSEnrichment()
function in Signac.
brain <- TSSEnrichment(brain, fast = FALSE)
brain$high.tss <- ifelse(brain$TSS.enrichment > 2, 'High', 'Low')
TSSPlot(brain, group.by = 'high.tss') + NoLegend()
brain$pct_reads_in_peaks <- brain$peak_region_fragments / brain$passed_filters * 100
brain$blacklist_ratio <- brain$blacklist_region_fragments / brain$peak_region_fragments
VlnPlot(
object = brain,
features = c('pct_reads_in_peaks', 'peak_region_fragments',
'TSS.enrichment', 'blacklist_ratio', 'nucleosome_signal'),
pt.size = 0.1,
ncol = 5
)
We remove cells that are outliers for these QC metrics.
brain <- subset(
x = brain,
subset = peak_region_fragments > 3000 &
peak_region_fragments < 100000 &
pct_reads_in_peaks > 40 &
blacklist_ratio < 0.025 &
nucleosome_signal < 4 &
TSS.enrichment > 2
)
brain
## An object of class Seurat
## 157203 features across 3512 samples within 1 assay
## Active assay: peaks (157203 features, 0 variable features)
brain <- RunTFIDF(brain)
brain <- FindTopFeatures(brain, min.cutoff = 'q0')
brain <- RunSVD(object = brain)
The first LSI component often captures sequencing depth (technical variation) rather than biological variation. If this is the case, the component should be removed from downstream analysis. We can assess the correlation between each LSI component and sequencing depth using the DepthCor()
function:
DepthCor(brain)
Here we see there is a very strong correlation between the first LSI component and the total number of counts for the cell, so we will perform downstream steps without this component.
Now that the cells are embedded in a low-dimensional space, we can use methods commonly applied for the analysis of scRNA-seq data to perform graph-based clustering, and non-linear dimension reduction for visualization. The functions RunUMAP()
, FindNeighbors()
, and FindClusters()
all come from the Seurat package.
brain <- RunUMAP(
object = brain,
reduction = 'lsi',
dims = 2:30
)
brain <- FindNeighbors(
object = brain,
reduction = 'lsi',
dims = 2:30
)
brain <- FindClusters(
object = brain,
algorithm = 3,
resolution = 1.2,
verbose = FALSE
)
DimPlot(object = brain, label = TRUE) + NoLegend()
# compute gene activities
gene.activities <- GeneActivity(brain)
# add the gene activity matrix to the Seurat object as a new assay
brain[['RNA']] <- CreateAssayObject(counts = gene.activities)
brain <- NormalizeData(
object = brain,
assay = 'RNA',
normalization.method = 'LogNormalize',
scale.factor = median(brain$nCount_RNA)
)
DefaultAssay(brain) <- 'RNA'
FeaturePlot(
object = brain,
features = c('Sst','Pvalb',"Gad2","Neurod6","Rorb","Syt6"),
pt.size = 0.1,
max.cutoff = 'q95',
ncol = 3
)
To help interpret the scATAC-seq data, we can classify cells based on an scRNA-seq experiment from the same biological system (the adult mouse brain). We utilize methods for cross-modality integration and label transfer, described here, with a more in-depth tutorial here.
You can download the raw data for this experiment from the Allen Institute website, and view the code used to construct this object on GitHub. Alternatively, you can download the pre-processed Seurat object here.
# Load the pre-processed scRNA-seq data
allen_rna <- readRDS("../vignette_data/allen_brain.rds")
allen_rna <- FindVariableFeatures(
object = allen_rna,
nfeatures = 5000
)
transfer.anchors <- FindTransferAnchors(
reference = allen_rna,
query = brain,
reduction = 'cca',
dims = 1:40
)
predicted.labels <- TransferData(
anchorset = transfer.anchors,
refdata = allen_rna$subclass,
weight.reduction = brain[['lsi']],
dims = 2:30
)
brain <- AddMetaData(object = brain, metadata = predicted.labels)
plot1 <- DimPlot(allen_rna, group.by = 'subclass', label = TRUE, repel = TRUE) + NoLegend() + ggtitle('scRNA-seq')
plot2 <- DimPlot(brain, group.by = 'predicted.id', label = TRUE, repel = TRUE) + NoLegend() + ggtitle('scATAC-seq')
plot1 + plot2
Why did we change default parameters?
We changed default parameters for FindIntegrationAnchors()
and FindVariableFeatures()
(including more features and dimensions). You can run the analysis both ways, and observe very similar results. However, when using default parameters we mislabel cluster 11 cells as Vip-interneurons, when they are in fact a Meis2 expressing CGE-derived interneuron population recently described by us and others. The reason is that this subset is exceptionally rare in the scRNA-seq data (0.3%), and so the genes define this subset (for example, Meis2) were too lowly expressed to be selected in the initial set of variable features. We therefore need more genes and dimensions to facilitate cross-modality mapping. Interestingly, this subset is 10-fold more abundant in the scATAC-seq data compared to the scRNA-seq data (see this paper for possible explanations.)
You can see that the RNA-based classifications are entirely consistent with the UMAP visualization, computed only on the ATAC-seq data. We can now easily annotate our scATAC-seq derived clusters (alternately, we could use the RNA classifications themselves). We note three small clusters (13, 20, 21) which represent subdivisions of the scRNA-seq labels. Try transferring the cluster label (which shows finer distinctions) from the allen scRNA-seq dataset, to annotate them!
Here, we find differentially accessible regions between excitatory neurons in different layers of the cortex.
#switch back to working with peaks instead of gene activities
DefaultAssay(brain) <- 'peaks'
da_peaks <- FindMarkers(
object = brain,
ident.1 = c("L2/3 IT"),
ident.2 = c("L4", "L5 IT", "L6 IT"),
min.pct = 0.05,
test.use = 'LR',
latent.vars = 'peak_region_fragments'
)
head(da_peaks)
## p_val avg_log2FC pct.1 pct.2 p_val_adj
## chr4-86523678-86525285 3.413350e-67 0.4663591 0.401 0.036 5.365888e-62
## chr15-87605281-87607659 2.945028e-59 0.4420617 0.485 0.090 4.629672e-54
## chr2-118700082-118704897 5.135950e-57 0.3178506 0.622 0.180 8.073868e-52
## chr4-101303935-101305131 7.414490e-54 0.4582150 0.340 0.025 1.165580e-48
## chr3-137056475-137058371 9.571201e-53 0.3479996 0.523 0.103 1.504621e-47
## chr13-69329933-69331707 1.019300e-52 -0.3604404 0.143 0.445 1.602370e-47
plot1 <- VlnPlot(
object = brain,
features = rownames(da_peaks)[1],
pt.size = 0.1,
idents = c("L4","L5 IT","L2/3 IT")
)
plot2 <- FeaturePlot(
object = brain,
features = rownames(da_peaks)[1],
pt.size = 0.1,
max.cutoff = 'q95'
)
plot1 | plot2
open_l23 <- rownames(da_peaks[da_peaks$avg_log2FC > 0.25, ])
open_l456 <- rownames(da_peaks[da_peaks$avg_log2FC < -0.25, ])
closest_l23 <- ClosestFeature(brain, open_l23)
closest_l456 <- ClosestFeature(brain, open_l456)
head(closest_l23)
## tx_id gene_name gene_id
## ENSMUST00000151481 ENSMUST00000151481 Fam154a ENSMUSG00000028492
## ENSMUSE00000647021 ENSMUST00000068088 Fam19a5 ENSMUSG00000054863
## ENSMUST00000104937 ENSMUST00000104937 Ankrd63 ENSMUSG00000078137
## ENSMUST00000131864 ENSMUST00000131864 Gm12796 ENSMUSG00000085721
## ENSMUST00000070198 ENSMUST00000070198 Ppp3ca ENSMUSG00000028161
## ENSMUST00000165341 ENSMUST00000165341 Otogl ENSMUSG00000091455
## gene_biotype type closest_region
## ENSMUST00000151481 protein_coding gap chr4-86487920-86538964
## ENSMUSE00000647021 protein_coding exon chr15-87625230-87625486
## ENSMUST00000104937 protein_coding cds chr2-118702266-118703438
## ENSMUST00000131864 lincRNA gap chr4-101292521-101318425
## ENSMUST00000070198 protein_coding utr chr3-136935226-136937727
## ENSMUST00000165341 protein_coding utr chr10-107762223-107762309
## query_region distance
## ENSMUST00000151481 chr4-86523678-86525285 0
## ENSMUSE00000647021 chr15-87605281-87607659 17570
## ENSMUST00000104937 chr2-118700082-118704897 0
## ENSMUST00000131864 chr4-101303935-101305131 0
## ENSMUST00000070198 chr3-137056475-137058371 118747
## ENSMUST00000165341 chr10-107751762-107753240 8982
head(closest_l456)
## tx_id gene_name gene_id
## ENSMUST00000044081 ENSMUST00000044081 Papd7 ENSMUSG00000034575
## ENSMUST00000042025 ENSMUST00000042025 Antxr1 ENSMUSG00000033420
## ENSMUST00000056102 ENSMUST00000056102 Dscam ENSMUSG00000050272
## ENSMUST00000189866 ENSMUST00000189866 Iqcf3 ENSMUSG00000023577
## ENSMUST00000163093 ENSMUST00000163093 Rnls ENSMUSG00000071573
## ENSMUST00000037205 ENSMUST00000037205 Mcee ENSMUSG00000033429
## gene_biotype type closest_region
## ENSMUST00000044081 protein_coding utr chr13-69497959-69499915
## ENSMUST00000042025 protein_coding gap chr6-87255882-87269242
## ENSMUST00000056102 protein_coding gap chr16-96629399-96640737
## ENSMUST00000189866 protein_coding gap chr9-106527019-106544398
## ENSMUST00000163093 protein_coding gap chr19-33202649-33382366
## ENSMUST00000037205 protein_coding utr chr7-64411988-64412121
## query_region distance
## ENSMUST00000044081 chr13-69329933-69331707 166251
## ENSMUST00000042025 chr6-87257016-87258706 0
## ENSMUST00000056102 chr16-96633303-96636474 0
## ENSMUST00000189866 chr9-106534540-106536125 0
## ENSMUST00000163093 chr19-33316841-33319335 0
## ENSMUST00000037205 chr7-64450755-64452725 38633
We can also create coverage plots grouped by cluster, cell type, or any other metadata stored in the object for any genomic region using the CoveragePlot()
function. These represent pseudo-bulk accessibility tracks, where signal from all cells within a group have been averaged together to visualize the DNA accessibility in a region.
# set plotting order
levels(brain) <- c("L2/3 IT","L4","L5 IT","L5 PT","L6 CT", "L6 IT","NP","Sst","Pvalb","Vip","Lamp5","Meis2","Oligo","Astro","Endo","VLMC","Macrophage")
CoveragePlot(
object = brain,
region = c("Neurod6", "Gad2"),
extend.upstream = 1000,
extend.downstream = 1000,
ncol = 1
)
Session Info
## R version 4.1.2 (2021-11-01)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 20.04.4 LTS
##
## Matrix products: default
## BLAS: /usr/lib/x86_64-linux-gnu/blas/libblas.so.3.9.0
## LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.9.0
##
## 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 stats graphics grDevices utils datasets methods
## [8] base
##
## other attached packages:
## [1] patchwork_1.1.1 ggplot2_3.3.5
## [3] EnsDb.Mmusculus.v79_2.99.0 ensembldb_2.18.2
## [5] AnnotationFilter_1.18.0 GenomicFeatures_1.46.1
## [7] AnnotationDbi_1.56.2 Biobase_2.54.0
## [9] GenomicRanges_1.46.1 GenomeInfoDb_1.30.0
## [11] IRanges_2.28.0 S4Vectors_0.32.2
## [13] BiocGenerics_0.40.0 SeuratObject_4.0.4
## [15] Seurat_4.1.0 Signac_1.6.0
##
## loaded via a namespace (and not attached):
## [1] rappdirs_0.3.3 SnowballC_0.7.0
## [3] rtracklayer_1.54.0 scattermore_0.8
## [5] ragg_1.2.0 tidyr_1.2.0
## [7] bit64_4.0.5 knitr_1.36
## [9] irlba_2.3.5 DelayedArray_0.20.0
## [11] data.table_1.14.2 rpart_4.1.16
## [13] KEGGREST_1.34.0 RCurl_1.98-1.5
## [15] generics_0.1.2 cowplot_1.1.1
## [17] RSQLite_2.2.8 RANN_2.6.1
## [19] future_1.24.0 bit_4.0.4
## [21] spatstat.data_2.1-2 xml2_1.3.2
## [23] httpuv_1.6.5 SummarizedExperiment_1.24.0
## [25] assertthat_0.2.1 xfun_0.28
## [27] hms_1.1.1 jquerylib_0.1.4
## [29] evaluate_0.14 promises_1.2.0.1
## [31] fansi_1.0.2 restfulr_0.0.13
## [33] progress_1.2.2 dbplyr_2.1.1
## [35] igraph_1.2.11 DBI_1.1.1
## [37] htmlwidgets_1.5.4 sparsesvd_0.2
## [39] spatstat.geom_2.3-2 purrr_0.3.4
## [41] ellipsis_0.3.2 RSpectra_0.16-0
## [43] dplyr_1.0.8 backports_1.3.0
## [45] biomaRt_2.50.0 deldir_1.0-6
## [47] MatrixGenerics_1.6.0 vctrs_0.3.8
## [49] Cairo_1.5-12.2 ROCR_1.0-11
## [51] abind_1.4-5 cachem_1.0.6
## [53] withr_2.5.0 ggforce_0.3.3
## [55] BSgenome_1.62.0 checkmate_2.0.0
## [57] sctransform_0.3.3 GenomicAlignments_1.30.0
## [59] prettyunits_1.1.1 goftest_1.2-3
## [61] cluster_2.1.2 lazyeval_0.2.2
## [63] crayon_1.5.0 hdf5r_1.3.5
## [65] pkgconfig_2.0.3 slam_0.1-50
## [67] labeling_0.4.2 tweenr_1.0.2
## [69] nlme_3.1-153 vipor_0.4.5
## [71] ProtGenerics_1.26.0 nnet_7.3-17
## [73] rlang_1.0.2 globals_0.14.0
## [75] lifecycle_1.0.1 miniUI_0.1.1.1
## [77] filelock_1.0.2 BiocFileCache_2.2.0
## [79] dichromat_2.0-0 ggrastr_1.0.0
## [81] rprojroot_2.0.2 polyclip_1.10-0
## [83] matrixStats_0.61.0 lmtest_0.9-39
## [85] Matrix_1.4-0 ggseqlogo_0.1
## [87] zoo_1.8-9 base64enc_0.1-3
## [89] beeswarm_0.4.0 ggridges_0.5.3
## [91] png_0.1-7 viridisLite_0.4.0
## [93] rjson_0.2.20 bitops_1.0-7
## [95] KernSmooth_2.23-20 Biostrings_2.62.0
## [97] blob_1.2.2 stringr_1.4.0
## [99] parallelly_1.30.0 spatstat.random_2.1-0
## [101] jpeg_0.1-9 scales_1.1.1
## [103] memoise_2.0.0 magrittr_2.0.2
## [105] plyr_1.8.6 ica_1.0-2
## [107] zlibbioc_1.40.0 compiler_4.1.2
## [109] BiocIO_1.4.0 RColorBrewer_1.1-2
## [111] fitdistrplus_1.1-6 Rsamtools_2.10.0
## [113] cli_3.2.0 XVector_0.34.0
## [115] listenv_0.8.0 pbapply_1.5-0
## [117] htmlTable_2.3.0 Formula_1.2-4
## [119] MASS_7.3-55 mgcv_1.8-39
## [121] tidyselect_1.1.2 stringi_1.7.6
## [123] textshaping_0.3.6 highr_0.9
## [125] yaml_2.3.5 latticeExtra_0.6-29
## [127] ggrepel_0.9.1 grid_4.1.2
## [129] sass_0.4.0 VariantAnnotation_1.40.0
## [131] fastmatch_1.1-3 tools_4.1.2
## [133] future.apply_1.8.1 parallel_4.1.2
## [135] rstudioapi_0.13 foreign_0.8-82
## [137] lsa_0.73.2 gridExtra_2.3
## [139] farver_2.1.0 Rtsne_0.15
## [141] digest_0.6.29 shiny_1.7.1
## [143] qlcMatrix_0.9.7 Rcpp_1.0.8
## [145] later_1.3.0 RcppAnnoy_0.0.19
## [147] httr_1.4.2 biovizBase_1.42.0
## [149] colorspace_2.0-3 XML_3.99-0.8
## [151] fs_1.5.2 tensor_1.5
## [153] reticulate_1.24 splines_4.1.2
## [155] uwot_0.1.11 RcppRoll_0.3.0
## [157] spatstat.utils_2.3-0 pkgdown_1.6.1.9001
## [159] plotly_4.10.0 systemfonts_1.0.3
## [161] xtable_1.8-4 jsonlite_1.8.0
## [163] R6_2.5.1 Hmisc_4.6-0
## [165] pillar_1.7.0 htmltools_0.5.2
## [167] mime_0.12 glue_1.6.2
## [169] fastmap_1.1.0 BiocParallel_1.28.1
## [171] codetools_0.2-18 utf8_1.2.2
## [173] lattice_0.20-45 bslib_0.3.1
## [175] spatstat.sparse_2.1-0 tibble_3.1.6
## [177] curl_4.3.2 ggbeeswarm_0.6.0
## [179] leiden_0.3.9 survival_3.2-13
## [181] rmarkdown_2.11 docopt_0.7.1
## [183] desc_1.4.0 munsell_0.5.0
## [185] GenomeInfoDbData_1.2.7 reshape2_1.4.4
## [187] gtable_0.3.0 spatstat.core_2.4-0