--- title: Analyzing single-cell RNA-seq data containing UMI counts author: - name: Aaron T. L. Lun affiliation: &CRUK Cancer Research UK Cambridge Institute, Li Ka Shing Centre, Robinson Way, Cambridge CB2 0RE, United Kingdom - name: Davis J. McCarthy affiliation: - &EMBL EMBL European Bioinformatics Institute, Wellcome Genome Campus, Hinxton, Cambridge CB10 1SD, United Kingdom - St Vincent's Institute of Medical Research, 41 Victoria Parade, Fitzroy, Victoria 3065, Australia - name: John C. Marioni affiliation: - *CRUK - *EMBL - Wellcome Trust Sanger Institute, Wellcome Genome Campus, Hinxton, Cambridge CB10 1SA, United Kingdom date: "`r Sys.Date()`" vignette: > %\VignetteIndexEntry{Analyzing scRNA-seq UMI count data} %\VignetteEngine{knitr::rmarkdown} %\VignetteEncoding{UTF-8} output: BiocStyle::html_document bibliography: ref.bib --- ```{r style, echo=FALSE, results='hide', message=FALSE} library(BiocStyle) library(knitr) opts_chunk$set(error=FALSE, message=FALSE, warning=FALSE) opts_chunk$set(fig.asp=1) ``` ```{r, echo=FALSE, results='hide'} all.urls <- c("https://storage.googleapis.com/linnarsson-lab-www-blobs/blobs/cortex/expression_mRNA_17-Aug-2014.txt", "https://storage.googleapis.com/linnarsson-lab-www-blobs/blobs/cortex/expression_mito_17-Aug-2014.txt", "https://storage.googleapis.com/linnarsson-lab-www-blobs/blobs/cortex/expression_spikes_17-Aug-2014.txt") all.basenames <- basename(all.urls) for (x in seq_along(all.urls)) { if (!file.exists(all.basenames[x])) { download.file(all.urls[x], all.basenames[x]) } } ``` # Overview In this workflow, we examine a heterogeneous dataset from a study of cell types in the mouse brain [@zeisel2015brain]. This contains approximately 3000 cells of varying types such as oligodendrocytes, microglia and neurons. Individual cells were isolated using the Fluidigm C1 microfluidics system [@pollen2014lowcoverage] and library preparation was performed on each cell using a UMI-based protocol. After sequencing, expression was quantified by counting the number of UMIs mapped to each gene. Count data for all endogenous genes, mitochondrial genes and spike-in transcripts were obtained from http://linnarssonlab.org/cortex. # Setting up the data The count data are distributed across several files, so some work is necessary to consolidate them into a single matrix. We define a simple utility function for loading data in from each file. (We stress that this function is only relevant to the current dataset, and should not be used for other datasets. This kind of effort is generally not required if all of the counts are in a single file and separated from the metadata.) ```{r} readFormat <- function(infile) { # First column is empty. metadata <- read.delim(infile, stringsAsFactors=FALSE, header=FALSE, nrow=10)[,-1] rownames(metadata) <- metadata[,1] metadata <- metadata[,-1] metadata <- as.data.frame(t(metadata)) # First column after row names is some useless filler. counts <- read.delim(infile, stringsAsFactors=FALSE, header=FALSE, row.names=1, skip=11)[,-1] counts <- as.matrix(counts) return(list(metadata=metadata, counts=counts)) } ``` Using this function, we read in the counts for the endogenous genes, ERCC spike-in transcripts and mitochondrial genes. ```{r} endo.data <- readFormat("expression_mRNA_17-Aug-2014.txt") spike.data <- readFormat("expression_spikes_17-Aug-2014.txt") mito.data <- readFormat("expression_mito_17-Aug-2014.txt") ``` We also need to rearrange the columns for the mitochondrial data, as the order is not consistent with the other files. ```{r} m <- match(endo.data$metadata$cell_id, mito.data$metadata$cell_id) mito.data$metadata <- mito.data$metadata[m,] mito.data$counts <- mito.data$counts[,m] ``` ```{r, echo=FALSE} stopifnot(identical(endo.data$metadata$cell_id, spike.data$metadata$cell_id)) # should be the same. stopifnot(all(endo.data$metadata$cell_id==mito.data$metadata$cell_id)) # should now be the same. ``` In this particular dataset, some genes are represented by multiple rows corresponding to alternative genomic locations. We sum the counts for all rows corresponding to a single gene for ease of interpretation. ```{r} raw.names <- sub("_loc[0-9]+$", "", rownames(endo.data$counts)) new.counts <- rowsum(endo.data$counts, group=raw.names, reorder=FALSE) endo.data$counts <- new.counts ``` The counts are then combined into a single matrix for constructing a `SingleCellExperiment` object. For convenience, metadata for all cells are stored in the same object for later access. ```{r} library(SingleCellExperiment) all.counts <- rbind(endo.data$counts, mito.data$counts, spike.data$counts) sce <- SingleCellExperiment(list(counts=all.counts), colData=endo.data$metadata) dim(sce) ``` We add gene-based annotation identifying rows that correspond to each class of features. We also determine the Ensembl identifier for each row. ```{r} # Specifying the nature of each row. nrows <- c(nrow(endo.data$counts), nrow(mito.data$counts), nrow(spike.data$counts)) is.spike <- rep(c(FALSE, FALSE, TRUE), nrows) is.mito <- rep(c(FALSE, TRUE, FALSE), nrows) isSpike(sce, "Spike") <- is.spike # Adding Ensembl IDs. library(org.Mm.eg.db) ensembl <- mapIds(org.Mm.eg.db, keys=rownames(sce), keytype="SYMBOL", column="ENSEMBL") rowData(sce)$ENSEMBL <- ensembl sce ``` ```{r, echo=FALSE, results='hide'} # Save some memory. rm(mito.data, endo.data, spike.data, new.counts) gc() ``` # Quality control on the cells The original authors of the study have already removed low-quality cells prior to data publication. Nonetheless, we compute some quality control metrics with `r Biocpkg("scater")` [@mccarthy2017scater] to check whether the remaining cells are satisfactory. ```{r} library(scater) sce <- calculateQCMetrics(sce, feature_controls=list(Mt=is.mito)) ``` We examine the distribution of the QC metrics across all cells (Figure \@ref(fig:libplotbrain)). The library sizes here are at least one order of magnitude lower than observed in the 416B dataset. This is consistent with the use of UMI counts rather than read counts, as each transcript molecule can only produce one UMI count but can yield many reads after fragmentation. In addition, the spike-in proportions are more variable than observed in the 416B dataset. This may reflect a greater variability in the total amount of endogenous RNA per cell when many cell types are present. ```{r libplotbrain, fig.wide=TRUE, fig.cap="Histograms of QC metrics including the library sizes, number of expressed genes and proportion of UMIs assigned to spike-in transcripts or mitochondrial genes for all cells in the brain dataset."} par(mfrow=c(2,2), mar=c(5.1, 4.1, 0.1, 0.1)) hist(sce$total_counts/1e3, xlab="Library sizes (thousands)", main="", breaks=20, col="grey80", ylab="Number of cells") hist(sce$total_features, xlab="Number of expressed genes", main="", breaks=20, col="grey80", ylab="Number of cells") hist(sce$pct_counts_Spike, xlab="ERCC proportion (%)", ylab="Number of cells", breaks=20, main="", col="grey80") hist(sce$pct_counts_Mt, xlab="Mitochondrial proportion (%)", ylab="Number of cells", breaks=20, main="", col="grey80") ``` We remove small outliers for the library size and the number of expressed features, and large outliers for the spike-in proportions. Again, the presence of spike-in transcripts means that we do not have to use the mitochondrial proportions. ```{r} libsize.drop <- isOutlier(sce$total_counts, nmads=3, type="lower", log=TRUE) feature.drop <- isOutlier(sce$total_features, nmads=3, type="lower", log=TRUE) spike.drop <- isOutlier(sce$pct_counts_Spike, nmads=3, type="higher") ``` Removal of low-quality cells is then performed by combining the filters for all of the metrics. The majority of cells are retained, which suggests that the original quality control procedures were generally adequate. ```{r} sce <- sce[,!(libsize.drop | feature.drop | spike.drop)] data.frame(ByLibSize=sum(libsize.drop), ByFeature=sum(feature.drop), BySpike=sum(spike.drop), Remaining=ncol(sce)) ``` We could improve our cell filtering procedure further by setting `batch` in `isOutlier` to one or more known factors, e.g., mouse/plate of origin. As previously mentioned, this would avoid inflation of the MAD and improve power to remove low-quality cells. However, for simplicity, we will not do this as sufficient quality control has already been performed. ```{r echo=FALSE, results='hide'} gc() ``` # Cell cycle classification Application of `cyclone` [@scialdone2015computational] to the brain dataset suggests that most of the cells are in G1 phase (Figure \@ref(fig:phaseplotbrain)). This requires the use of the Ensembl identifiers to match up with the pre-defined classifier. ```{r phaseplotbrain, message=FALSE, fig.cap="Cell cycle phase scores from applying the pair-based classifier on the brain dataset, where each point represents a cell."} library(scran) mm.pairs <- readRDS(system.file("exdata", "mouse_cycle_markers.rds", package="scran")) assignments <- cyclone(sce, mm.pairs, gene.names=rowData(sce)$ENSEMBL) table(assignments$phase) plot(assignments$score$G1, assignments$score$G2M, xlab="G1 score", ylab="G2/M score", pch=16) ``` However, the intepretation of this result requires some caution due to differences between the training and test datasets. The classifier was trained on C1 SMARTer data and accounts for the biases in that protocol. The brain dataset uses UMI counts, which has a different set of biases, e.g., 3'-end coverage only, no length bias, no amplification noise. Furthermore, many neuronal cell types are expected to lie in the G0 resting phase, which is distinct from the other phases of the cell cycle [@coller2006new]. `cyclone` will generally assign such cells to the closest known phase in the training set, which would be G1. ```{r echo=FALSE, results='hide'} gc() ``` # Examining gene-level metrics Figure \@ref(fig:topgenebrain) shows the most highly expressed genes across the cell population in the brain dataset. This is mostly occupied by spike-in transcripts, reflecting the use of spike-in concentrations that span the entire range of expression. There are also a number of constitutively expressed genes, as expected. ```{r topgenebrain, fig.asp=1.2, fig.wide=TRUE, fig.cap="Percentage of total counts assigned to the top 50 most highly-abundant features in the brain dataset. For each feature, each bar represents the percentage assigned to that feature for a single cell, while the circle represents the average across all cells. Bars are coloured by the total number of expressed features in each cell, while circles are coloured according to whether the feature is labelled as a control feature."} fontsize <- theme(axis.text=element_text(size=12), axis.title=element_text(size=16)) plotQC(sce, type = "highest-expression", n=50) + fontsize ``` Gene abundance is quantified by computing the average count across all cells (Figure \@ref(fig:abhistbrain)). As previously mentioned, the UMI count is generally lower than the read count. ```{r abhistbrain, fig.cap="Histogram of log-average counts for all genes in the brain dataset."} ave.counts <- calcAverage(sce, use_size_factors=FALSE) hist(log10(ave.counts), breaks=100, main="", col="grey", xlab=expression(Log[10]~"average count")) ``` We save the average counts into the `SingleCellExperiment` object for later use. We also remove genes that have average counts of zero, as this means that they are not expressed in any cell. ```{r} rowData(sce)$ave.count <- ave.counts to.keep <- ave.counts > 0 sce <- sce[to.keep,] summary(to.keep) ``` ```{r echo=FALSE, results='hide'} gc() ``` # Normalization of cell-specific biases For endogenous genes, normalization is performed using the `computeSumFactors` function as previously described. Here, we cluster similar cells together and normalize the cells in each cluster using the deconvolution method [@lun2016pooling]. This improves normalization accuracy by reducing the number of DE genes between cells in the same cluster. Scaling is then performed to ensure that size factors of cells in different clusters are comparable. ```{r} clusters <- quickCluster(sce, min.mean=0.1, method="igraph") sce <- computeSumFactors(sce, cluster=clusters, min.mean=0.1) summary(sizeFactors(sce)) ``` We use a average count threshold of 0.1 to define high-abundance genes to use during normalization. This is lower than the default threshold of `min.mean=1` in `computeSumFactors`, reflecting the fact that UMI counts are generally smaller than read counts. ```{r echo=FALSE, results='hide'} gc() ``` Compared to the 416B analysis, more scatter is observed around the trend between the total count and size factor for each cell (Figure \@ref(fig:normplotbrain)). This is consistent with an increased amount of DE between cells of different types, which compromises the accuracy of library size normalization [@robinson2010scaling]. In contrast, the size factors are estimated based on median ratios and are more robust to the presence of DE between cells. ```{r normplotbrain, fig.cap="Size factors from deconvolution, plotted against library sizes for all cells in the brain dataset. Axes are shown on a log-scale."} plot(sizeFactors(sce), sce$total_counts/1e3, log="xy", ylab="Library size (thousands)", xlab="Size factor") ``` We also compute size factors specific to the spike-in set, as previously described. ```{r} sce <- computeSpikeFactors(sce, type="Spike", general.use=FALSE) ``` Finally, normalized log-expression values are computed for each endogenous gene or spike-in transcript using the appropriate size factors. ```{r} sce <- normalize(sce) ``` ```{r echo=FALSE, results='hide'} gc() ``` __Comments from Aaron:__ - Only a rough clustering is required to avoid pooling together very different cell types in `computeSumFactors`. The function is robust to a moderate level of differential expression between cells in the same cluster. - `quickCluster` uses distances based on Spearman's rank correlation for clustering. This ensures that scaling biases in the counts do not affect clustering, but yields very coarse clusters and is not recommended for biological interpretation. - For large datasets, using `method="igraph"` in `quickCluster` will speed up clustering. This uses a graph-based clustering algorithm - see `?buildSNNGraph` for more details. # Modelling and removing technical noise We model the technical noise by fitting a mean-variance trend to the spike-in transcripts, as previously described. In theory, we should block on the plate of origin for each cell. However, only 20-40 cells are available on each plate, and the population is also highly heterogeneous. This means that we cannot assume that the distribution of sampled cell types on each plate is the same. Thus, to avoid regressing out potential biology, we will not block on any factors in this analysis. ```{r} var.fit <- trendVar(sce, parametric=TRUE, loess.args=list(span=0.4)) var.out <- decomposeVar(sce, var.fit) ``` Figure \@ref(fig:hvgplotbrain) indicates that the trend is fitted accurately to the technical variances. The technical and total variances are also much smaller than those in the 416B dataset. This is due to the use of UMIs, which reduces the noise caused by variable PCR amplification [@islam2014quantitative]. Furthermore, the spike-in trend is consistently lower than the variances of the endogenous genes. This reflects the heterogeneity in gene expression across cells of different types. ```{r hvgplotbrain, fig.cap="Variance of normalized log-expression values against the mean for each gene, calculated across all cells in the brain dataset after blocking on the sex effect. The blue line represents the mean-dependent trend in the technical variance of the spike-in transcripts (also highlighted as red points)."} plot(var.out$mean, var.out$total, pch=16, cex=0.6, xlab="Mean log-expression", ylab="Variance of log-expression") points(var.out$mean[isSpike(sce)], var.out$total[isSpike(sce)], col="red", pch=16) curve(var.fit$trend(x), col="dodgerblue", add=TRUE, lwd=2) ``` We check the distribution of expression values for the genes with the largest biological components to ensure that they are not driven by outliers (Figure \@ref(fig:hvgvioplotbrain)). Some tweaking of the `plotExpression` parameters is necessary to visualize a large number of cells. ```{r hvgvioplotbrain, fig.cap="Violin plots of normalized log-expression values for the top 10 HVGs in the brain dataset. For each gene, each point represents the log-expression value for an individual cell."} chosen.genes <- order(var.out$bio, decreasing=TRUE)[1:10] plotExpression(sce, rownames(var.out)[chosen.genes], alpha=0.05, jitter="jitter") + fontsize ``` Finally, we use PCA to denoise the expression values, yielding a set of coordinates for each cell where the technical noise has been removed. Setting `approximate=TRUE` in `denoisePCA` will perform an approximate singular value decomposition, using methods from the `r CRANpkg("irlba")` package. This is much faster than the exact algorithm on large datasets without much loss of accuracy. ```{r} sce <- denoisePCA(sce, technical=var.fit$trend, approximate=TRUE) ncol(reducedDim(sce, "PCA")) ``` ```{r, echo=FALSE, results='hide', message=FALSE} gc() ``` # Data exploration with dimensionality reduction We perform dimensionality reduction on the denoised PCs to check if there is any substructure. Cells separate into clear clusters in the _t_-SNE plot [@van2008visualizing] in Figure \@ref(fig:tsneplotbrain), corresponding to distinct subpopulations. This is consistent with the presence of multiple cell types in the diverse brain population. We increase the perplexity to favour visualization of the overall structure at the expense of local scale. ```{r tsneplotbrain, fig.cap="_t_-SNE plots constructed from the denoised PCs of the brain dataset. Each point represents a cell and is coloured according to its expression of _Neurod6_ (left) or _Mog_ (right).", fig.width=12, fig.asp=0.5} sce <- runTSNE(sce, use_dimred="PCA", perplexity=50, rand_seed=1000) tsne1 <- plotTSNE(sce, colour_by="Neurod6") + fontsize tsne2 <- plotTSNE(sce, colour_by="Mog") + fontsize multiplot(tsne1, tsne2, cols=2) ``` The PCA plot is less effective at separating cells into many different clusters (Figure \@ref(fig:pcaplotbrain)). This is because the first two PCs are driven by strong differences between specific subpopulations, which reduces the resolution of more subtle differences between some of the other subpopulations. Nonetheless, some substructure is still visible. ```{r pcaplotbrain, fig.cap="PCA plots constructed from the denoised PCs of the brain dataset. Each point represents a cell and is coloured according to its expression of the _Neurod6_ (left) or _Mog_ (right).", fig.width=12, fig.asp=0.5} pca1 <- plotReducedDim(sce, use_dimred="PCA", colour_by="Neurod6") + fontsize pca2 <- plotReducedDim(sce, use_dimred="PCA", colour_by="Mog") + fontsize multiplot(pca1, pca2, cols=2) ``` For both methods, we colour each cell based on the expression of a particular gene. This is a useful strategy for visualizing changes in expression across the lower-dimensional space. It can also be used to characterise each cluster if the selected genes are known markers for particular cell types. For example, _Mog_ can be used to identify clusters corresponding to oligodendrocytes. ```{r echo=FALSE, results='hide'} gc() ``` # Clustering cells into putative subpopulations The reduced dimension coordinates are used to cluster cells into putative subpopulations. We do so by constructing a shared-nearest-neighbour graph [@xu2015identification], in which cells are the nodes and edges are formed between cells that share nearest neighbours. Clusters are then defined as highly connected communities of cells within this graph, using methods from the `r CRANpkg("igraph")` package. This is more efficient than forming a pairwise distance matrix for hierarchical clustering of large numbers of cells. ```{r} snn.gr <- buildSNNGraph(sce, use.dimred="PCA") cluster.out <- igraph::cluster_walktrap(snn.gr) my.clusters <- cluster.out$membership table(my.clusters) ``` The modularity score provides a global measure of clustering performance for community detection methods. Briefly, it compares the number of within-cluster edges to the expected number under a null model of random edges. A high modularity score (approaching the maximum of 1) indicates that the detected clusters are enriched for internal edges, with relatively few edges between clusters. ```{r} igraph::modularity(cluster.out) ``` We further investigate the clusters by examining the total weight of edges for each pair of clusters. For each pair, the observed total weight is compared to what is expected under a null model, similar to the modularity calculation. Most clusters contain more internal links than expected (Figure \@ref(fig:heatmodbrain)), while links between clusters are fewer than expected. This indicates that we successfully clustered cells into highly-connected communities. ```{r heatmodbrain, fig.cap="Heatmap of the log~10~-ratio of the total weight between nodes in the same cluster or in different clusters, relative to the total weight expected under a null model of random links."} mod.out <- clusterModularity(snn.gr, my.clusters, get.values=TRUE) ratio <- log10(mod.out$observed/mod.out$expected + 1) library(pheatmap) pheatmap(ratio, cluster_rows=FALSE, cluster_cols=FALSE, color=colorRampPalette(c("white", "blue"))(100)) ``` We visualize the cluster assignments for all cells on the _t_-SNE plot in Figure \@ref(fig:tsneclusterbrain). Adjacent cells are generally assigned to the same cluster, indicating that the clustering procedure was applied correctly. ```{r tsneclusterbrain, message=FALSE, fig.cap="_t_-SNE plot of the denoised PCs of the brain dataset. Each point represents a cell and is coloured according to its assigned cluster identity."} sce$cluster <- factor(my.clusters) plotTSNE(sce, colour_by="cluster") + fontsize ``` ```{r echo=FALSE, results='hide'} gc() ``` An alternative approach is to use graph-based visualizations such as force-directed layouts (Figure \@ref(fig:fdlbrain)). These are appealing as they directly represent the relationships used during clustering. However, convergence tends to be slow for large graphs, so some tinkering with `niter=` may be required to ensure that the results are stable. ```{r fdlbrain, message=FALSE, fig.cap="Force-directed layout for the shared nearest-neighbour graph of the brain dataset. Each point represents a cell and is coloured according to its assigned cluster identity."} set.seed(2000) reducedDim(sce, "force") <- igraph::layout_with_fr(snn.gr, niter=5000) plotReducedDim(sce, colour_by="cluster", use_dimred="force") ``` Very heterogeneous datasets may yield a few large clusters on the first round of clustering. It can be useful to repeat the variance modelling, denoising and clustering using only the cells within each of the initial clusters. This can be achieved by subsetting `sce` according to a particular level of `my.clusters`, and re-applying the relevant functions on the subset. Doing so may focus on a different set of genes that define heterogeneity _within_ an initial cluster, as opposed to those that define differences _between_ the initial clusters. This would allow fine-scale structure within each cluster to be explored at greater resolution. For simplicity, though, we will only use the broad clusters corresponding to clear subpopulations in this workflow. **Comments from Aaron:** - Many different clustering methods are available in the `r CRANpkg("igraph")` package. We find that the Walktrap algorithm is usually a good default choice [@yang2016comparative], though users are encouraged to experiment with different algorithms. - Decreasing the number of neighbours `k` in `buildSNNGraph` will reduce the connectivity of the graph. This will generally result in the formation of smaller clusters [@xu2015identification], which may be desirable if greater resolution is required. - We do not use the silhouette coefficient to assess clustering for large datasets. This is because `cluster::silhouette` requires the construction of a distance matrix, which may not be feasible when many cells are involved. - Notice that we do not run `library(igraph)`, but instead use `igraph::` to extract methods from the package. This is because `r CRANpkg("igraph")` contains a `normalize` method that will override its counterpart from `r Biocpkg("scater")`, resulting in some unusual bugs. # Detecting marker genes between subpopulations We use the `findMarkers` function with `direction="up"` to identify upregulated marker genes for each cluster. As previously mentioned, we focus on upregulated genes as these can quickly provide positive identification of cell type in a heterogeneous population. We examine the table for cluster 1, in which log-fold changes are reported between cluster 1 and every other cluster. The same output is provided for each cluster in order to identify genes that discriminate between clusters. ```{r, echo=FALSE, results="hide"} old.digits <- options()$digits options(digits=3) ``` ```{r} markers <- findMarkers(sce, my.clusters, direction="up") marker.set <- markers[["1"]] head(marker.set[,1:8], 10) # only first 8 columns, for brevity ``` ```{r, echo=FALSE, results="hide"} options(digits=old.digits) ``` We save the list of candidate marker genes for further examination, using compression to reduce the file size. The `overlapExprs` function may also be useful for summarizing differences between clusters, as previously mentioned. ```{r} gzout <- gzfile("brain_marker_1.tsv.gz", open="wb") write.table(marker.set, file=gzout, sep="\t", quote=FALSE, row.names=FALSE) close(gzout) ``` Figure \@ref(fig:heatmapmarkerbrain) indicates that most of the top markers are strongly DE in cells of cluster 1 compared to some or all of the other clusters. We can use these markers to identify cells from cluster 1 in validation studies with an independent population of cells. A quick look at the markers suggest that cluster 1 represents interneurons based on expression of *Gad1* and *Slc6a1* [@zeng2012largescale], differing from closely related cells in cluster 10 by virtue of high *Synpr* expression. ```{r heatmapmarkerbrain, fig.wide=TRUE, fig.cap="Heatmap of mean-centred and normalized log-expression values for the top set of markers for cluster 1 in the brain dataset. Column colours represent the cluster to which each cell is assigned, as indicated by the legend."} top.markers <- rownames(marker.set)[marker.set$Top <= 10] plotHeatmap(sce, features=top.markers, columns=order(my.clusters), colour_columns_by="cluster", cluster_cols=FALSE, center=TRUE, symmetric=TRUE, zlim=c(-5, 5)) ``` # Concluding remarks Having completed the basic analysis, we save the `SingleCellExperiment` object with its associated data to file. This is especially important here as the brain dataset is quite large. If further analyses are to be performed, it would be inconvenient to have to repeat all of the pre-processing steps described above. ```{r} saveRDS(file="brain_data.rds", sce) ``` ```{r, echo=FALSE, results='hide'} gc() ``` All software packages used in this workflow are publicly available from the Comprehensive R Archive Network (https://cran.r-project.org) or the Bioconductor project (http://bioconductor.org). The specific version numbers of the packages used are shown below, along with the version of the R installation. ```{r} sessionInfo() ``` # References