--- title: "Accessing the Human Cell Atlas datasets" author: - name: Federico Marini affiliation: - &id1 Institute of Medical Biostatistics, Epidemiology and Informatics (IMBEI), Mainz - Center for Thrombosis and Hemostasis (CTH), Mainz email: marinif@uni-mainz.de date: "`r BiocStyle::doc_date()`" package: "`r BiocStyle::pkg_ver('HCAData')`" output: BiocStyle::html_document: toc_float: true vignette: > %\VignetteIndexEntry{Accessing the Human Cell Atlas datasets} %\VignetteEngine{knitr::rmarkdown} %\VignetteEncoding{UTF-8} %\VignettePackage{HCAData} --- **Compiled date**: `r Sys.Date()` **Last edited**: 2019-05-24 **License**: `r packageDescription("HCAData")[["License"]]` ```{r setup, include = FALSE} knitr::opts_chunk$set( collapse = TRUE, comment = "#>", error = FALSE, warning = FALSE, message = FALSE ) ``` # Introducing the `HCAData` package The `r Biocpkg("HCAData")` package allows a direct access to the dataset generated by the Human Cell Atlas project for further processing in R and Bioconductor. It does so by providing the datasets as `r Biocpkg("SingleCellExperiment")` objects, i.e. a format which is both efficient and very widely adopted throughout many existing Bioconductor workflows. The datasets are otherwise available in other formats (also as raw data) at this link: http://preview.data.humancellatlas.org/. # Installing `HCAData` The `r Biocpkg("HCAData")` package can be installed in the conventional way via `r CRANpkg("BiocManager")`. ```{r install, eval=FALSE} if (!requireNamespace("BiocManager", quietly = TRUE)) install.packages("BiocManager") BiocManager::install("HCAData") ``` This package makes extensive use of the `r Biocpkg("HDF5Array")` package to avoid loading the entire data set in memory. Instead, it stores the counts on disk as a HDF5 file, and loads subsets of the data into memory upon request. # Loading `HCAData` and the Human Cell Atlas data ```{r loadpkg} library("HCAData") ``` We use the `HCAData` function to download the relevant files from Bioconductor's `r Biocpkg("ExperimentHub")` web resource. If no argument is provided, a list of the available datasets is returned, specifying which name to enter as `dataset` parameter when calling `HCAData`. ```{r firstcall} HCAData() ``` The list of relevant files includes the HDF5 file containing the counts, as well as the metadata on the rows (genes) and columns (cells). The output is a single `SingleCellExperiment` object from the `r Biocpkg("SingleCellExperiment")` package. Being based on `r Biocpkg("ExperimentHub")`, the data related to this package can be accessed and queried directly using the package name. Retrieval is then as easy as using their ExperimentHub accession numbers (for the single components of each set), or by using the convenience function provided in this package. ```{r fetchingdata} suppressPackageStartupMessages({ library(ExperimentHub) library(SingleCellExperiment) }) eh <- ExperimentHub() query(eh, "HCAData") # these three are the components to the bone marrow dataset bonemarrow_h5densematrix <- eh[["EH2047"]] bonemarrow_coldata <- eh[["EH2048"]] bonemarrow_rowdata <- eh[["EH2049"]] # and are put together when calling... sce_bonemarrow <- HCAData("ica_bone_marrow") sce_bonemarrow # similarly, to access the umbilical cord blood dataset sce_cordblood <- HCAData("ica_cord_blood") sce_cordblood ``` # Explore data with `iSEE` The datasets are provided in the form of a `SingleCellExperiment` object. A natural companion to this data structure is the `r Biocpkg("iSEE")` package, which can be used for interactive and reproducible data exploration. Any analysis steps should be performed in advance before calling `iSEE`, and since these datasets can be quite big, the operations can be time consuming, and/or require a considerable amount of resources. ## Processing a subset of the HCA bone marrow data For the scope of the vignette, we subset some cells in the bone marrow dataset to reduce the runtime, and apply some of the steps one would ideally follow in analysing droplet based datasets. For more information on how to properly process such datasets, please refer to the amazing set of resources available in the `r Biocpkg("simpleSingleCell")` workflow package. In brief: we start loading the required libraries for preprocessing, and taking a subset of the bone marrow dataset ```{r subset} library(scran) library(BiocSingular) library(scater) set.seed(42) sce <- sce_bonemarrow[, sample(seq_len(ncol(sce_bonemarrow)), 1000, replace = FALSE)] ``` First, we relabel the rows with the gene symbols for easier reading with the `uniquifyFeatureNames()` function from `r Biocpkg("scater")`. Then we compute some QC metrics and add these to the original `sce` object, using `addPerCellQC()`. ```{r featureanno} rownames(sce) <- uniquifyFeatureNames(rowData(sce)$ID, rowData(sce)$Symbol) head(rownames(sce)) counts(sce) <- as.matrix(counts(sce)) sce <- scater::addPerCellQC(sce) ``` We proceed with normalization, performed with the deconvolution method implemented in `r Biocpkg("scran")` - a pre-clustering step is done in advance, to avoid pooling cells that are very different between each other. ```{r clustnorm} lib.sf.bonemarrow <- librarySizeFactors(sce) summary(lib.sf.bonemarrow) set.seed(42) clusters <- quickCluster(sce) table(clusters) sce <- computeSumFactors(sce, min.mean=0.1, cluster=clusters) sce <- logNormCounts(sce) assayNames(sce) ``` In the following lines of code, we model the mean-variance trend (with technical noise as Poisson), to extract the biological component of the variance for the genes under inspection. ```{r hvg} dec.bonemarrow <- modelGeneVarByPoisson(sce) top.dec <- dec.bonemarrow[order(dec.bonemarrow$bio, decreasing=TRUE),] head(top.dec) fit.bonemarrow <- metadata(dec.bonemarrow) plot(fit.bonemarrow$mean, fit.bonemarrow$var, xlab="Mean of log-expression", ylab="Variance of log-expression", pch = 16) curve(fit.bonemarrow$trend(x), col="dodgerblue", add=TRUE, lwd=2) plotExpression(sce, features=rownames(top.dec)[1:10]) ``` With the `denoisePCA` function from `r Biocpkg("scran")`, we can compute a reduced dimension view of the data, accounting for the Poisson technical trend. Once we computed the PCA, we provide that as initialization to `runTSNE`, and obtain the t-SNE results, stored in the appropriate slot of our `sce` object. ```{r dimred} hvg.bonemarrow <- getTopHVGs(dec.bonemarrow, prop = 0.1) set.seed(42) sce <- denoisePCA(sce, technical=dec.bonemarrow, subset.row = hvg.bonemarrow, BSPARAM=IrlbaParam()) ncol(reducedDim(sce, "PCA")) plot(attr(reducedDim(sce), "percentVar"), xlab="PC", ylab="Proportion of variance explained") abline(v=ncol(reducedDim(sce, "PCA")), lty=2, col="red") plotPCA(sce, ncomponents=3, colour_by="percent_top_500") set.seed(42) sce <- runTSNE(sce, dimred="PCA", perplexity=30) plotTSNE(sce, colour_by="percent_top_500") ``` After this, we can compute a shared nearest neighbour graph to identify clusters in our dataset. Once the cluster memberships are defined, we assign this vector to a corresponding `colData` slot, and plot the t-SNE for our subset, coloured accordingly. ```{r clusters} snn.gr <- buildSNNGraph(sce, use.dimred="PCA") clusters <- igraph::cluster_walktrap(snn.gr) sce$Cluster <- factor(clusters$membership) table(sce$Cluster) plotTSNE(sce, colour_by="Cluster") ``` Even if we only took a subset of the available full data, we can observe that different clusters are indeed nicely separated. Subsequent steps would then involve identification of marker genes, and more advanced downstream techniques, which are not part of the scope of this vignette. To read more on what can be done, the vignettes of the `r Biocpkg("simpleSingleCell")` workflow package are an excellent place to start. ## Exploring the dataset with `iSEE` Once the processing steps above are done, we can call `iSEE` with the subsampled `SingleCellExperiment` object. ```{r launchisee, eval=FALSE} if (require(iSEE)) { iSEE(sce) } ``` ## Saving the processed object You can save the `sce` object to a serialized R object with ```{r save, eval=FALSE} destination <- "where/to/store/the/processed/data.rds" saveRDS(sce, file = destination) ``` The object can be read into a new R session with `readRDS(destination)`, provided the HDF5 file remains in its original location (conveniently stored in the default location of `ExperimentHub`). # Session info {-} ```{r sessioninfo} sessionInfo() ```