The Barcode, UMI, Set (BUS) format is a new way to represent pseudoalignments of reads from RNA-seq. Files of this format can be efficiently generated by the command line tool kallisto bus
. With kallisto bus
and this package, we go from the fastq files to the sparse matrix used for downstream analysis such as with Seurat within half an hour, while Cell Ranger would take hours.
In this vignette, we convert an 10x 1:1 mouse and human cell mixture dataset from the BUS format to a sparse matrix. To see how the BUS format can be generated from fastq file, as well as more in depth vignettes, see the website of this package.
Note that this vignette is deprecated and is kept for historical reasons as it was implemented when kallisto | bustools
was experimental. The functionality of make_sparse_matrix
has been implemented more efficiently in the command line tool bustools
. Please use the updated version of bustools
and if you wish, the wrapper kb
instead.
# The dataset package
library(TENxBUSData)
library(BUSpaRse)
library(Matrix)
library(zeallot)
library(ggplot2)
TENxBUSData(".", dataset = "hgmm100")
#> snapshotDate(): 2021-10-18
#> see ?TENxBUSData and browseVignettes('TENxBUSData') for documentation
#> loading from cache
#> The downloaded files are in /tmp/Rtmp605ceL/Rbuild35ba1e2bb78f6f/BUSpaRse/vignettes/out_hgmm100
#> [1] "/tmp/Rtmp605ceL/Rbuild35ba1e2bb78f6f/BUSpaRse/vignettes/out_hgmm100"
First, we map transcripts, as in the kallisto
index, to the corresponding genes.
tr2g <- transcript2gene(species = c("Homo sapiens", "Mus musculus"), type = "vertebrate",
kallisto_out_path = "./out_hgmm100", ensembl_version = 99,
write_tr2g = FALSE)
#> Querying biomart for transcript and gene IDs of Homo sapiens
#> Querying biomart for transcript and gene IDs of Mus musculus
head(tr2g)
#> transcript gene gene_name chromosome_name
#> 1 ENST00000434970.2 ENSG00000237235.2 TRDD2 14
#> 2 ENST00000415118.1 ENSG00000223997.1 TRDD1 14
#> 3 ENST00000448914.1 ENSG00000228985.1 TRDD3 14
#> 4 ENST00000632684.1 ENSG00000282431.1 TRBD1 7
#> 5 ENST00000390583.1 ENSG00000211923.1 IGHD3-10 14
#> 6 ENST00000431440.2 ENSG00000232543.2 IGHD4-11 14
Here we make both the gene count matrix and the TCC matrix.
c(gene_count, tcc) %<-% make_sparse_matrix("./out_hgmm100/output.sorted.txt",
tr2g = tr2g, est_ncells = 1e5,
est_ngenes = nrow(tr2g))
#> Reading matrix.ec
#> Processing ECs
#> Matching genes to ECs
#> Reading data
#> Constructing gene count matrix
#> Constructing TCC matrix
Here we have a sparse matrix with genes in rows and cells in columns.
dim(gene_count)
#> [1] 26937 151449
This dataset should only have about 100 cells, but here we get over 100,000. In fact, most of the barcodes correspond to empty droplets; they can be removed by filtering out barcodes with too few UMI.
tot_counts <- Matrix::colSums(gene_count)
summary(tot_counts)
#> Min. 1st Qu. Median Mean 3rd Qu. Max.
#> 1.00 1.00 2.00 24.77 5.00 64041.00
df1 <- get_knee_df(gene_count)
infl1 <- get_inflection(df1)
knee_plot(df1, infl1)
gene_count <- gene_count[, tot_counts > infl1]
dim(gene_count)
#> [1] 26937 122
Then this sparse matrix can be used in Seurat for downstream analysis.
Likewise, we can remove empty droplets from the TCC matrix.
dim(tcc)
#> [1] 129371 157659
This dataset should only have about 100 cells, but here we get over 100,000. In fact, most of the barcodes correspond to empty droplets; they can be removed by filtering out barcodes with too few UMI.
tot_counts <- Matrix::colSums(tcc)
summary(tot_counts)
#> Min. 1st Qu. Median Mean 3rd Qu. Max.
#> 1.00 1.00 2.00 25.84 5.00 69235.00
df2 <- get_knee_df(tcc)
infl2 <- get_inflection(df2)
knee_plot(df2, infl2)
tcc <- tcc[, tot_counts > infl2]
dim(tcc)
#> [1] 129371 121
sessionInfo()
#> R version 4.1.1 (2021-08-10)
#> Platform: x86_64-pc-linux-gnu (64-bit)
#> Running under: Ubuntu 20.04.3 LTS
#>
#> Matrix products: default
#> BLAS: /home/biocbuild/bbs-3.14-bioc/R/lib/libRblas.so
#> LAPACK: /home/biocbuild/bbs-3.14-bioc/R/lib/libRlapack.so
#>
#> 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
#>
#> attached base packages:
#> [1] stats graphics grDevices utils datasets methods base
#>
#> other attached packages:
#> [1] ggplot2_3.3.5 zeallot_0.1.0 Matrix_1.3-4 BUSpaRse_1.8.0
#> [5] TENxBUSData_1.7.0 BiocStyle_2.22.0
#>
#> loaded via a namespace (and not attached):
#> [1] ProtGenerics_1.26.0 bitops_1.0-7
#> [3] matrixStats_0.61.0 bit64_4.0.5
#> [5] progress_1.2.2 filelock_1.0.2
#> [7] httr_1.4.2 GenomeInfoDb_1.30.0
#> [9] tools_4.1.1 bslib_0.3.1
#> [11] utf8_1.2.2 R6_2.5.1
#> [13] colorspace_2.0-2 DBI_1.1.1
#> [15] lazyeval_0.2.2 BiocGenerics_0.40.0
#> [17] withr_2.4.2 prettyunits_1.1.1
#> [19] tidyselect_1.1.1 bit_4.0.4
#> [21] curl_4.3.2 compiler_4.1.1
#> [23] Biobase_2.54.0 xml2_1.3.2
#> [25] DelayedArray_0.20.0 rtracklayer_1.54.0
#> [27] bookdown_0.24 sass_0.4.0
#> [29] scales_1.1.1 rappdirs_0.3.3
#> [31] stringr_1.4.0 digest_0.6.28
#> [33] Rsamtools_2.10.0 rmarkdown_2.11
#> [35] XVector_0.34.0 pkgconfig_2.0.3
#> [37] htmltools_0.5.2 MatrixGenerics_1.6.0
#> [39] highr_0.9 ensembldb_2.18.0
#> [41] dbplyr_2.1.1 fastmap_1.1.0
#> [43] BSgenome_1.62.0 rlang_0.4.12
#> [45] RSQLite_2.2.8 shiny_1.7.1
#> [47] farver_2.1.0 jquerylib_0.1.4
#> [49] BiocIO_1.4.0 generics_0.1.1
#> [51] jsonlite_1.7.2 BiocParallel_1.28.0
#> [53] dplyr_1.0.7 RCurl_1.98-1.5
#> [55] magrittr_2.0.1 GenomeInfoDbData_1.2.7
#> [57] munsell_0.5.0 Rcpp_1.0.7
#> [59] S4Vectors_0.32.0 fansi_0.5.0
#> [61] lifecycle_1.0.1 stringi_1.7.5
#> [63] yaml_2.2.1 SummarizedExperiment_1.24.0
#> [65] zlibbioc_1.40.0 BiocFileCache_2.2.0
#> [67] AnnotationHub_3.2.0 grid_4.1.1
#> [69] blob_1.2.2 parallel_4.1.1
#> [71] promises_1.2.0.1 ExperimentHub_2.2.0
#> [73] crayon_1.4.1 lattice_0.20-45
#> [75] Biostrings_2.62.0 GenomicFeatures_1.46.0
#> [77] hms_1.1.1 KEGGREST_1.34.0
#> [79] magick_2.7.3 knitr_1.36
#> [81] pillar_1.6.4 GenomicRanges_1.46.0
#> [83] rjson_0.2.20 biomaRt_2.50.0
#> [85] stats4_4.1.1 XML_3.99-0.8
#> [87] glue_1.4.2 BiocVersion_3.14.0
#> [89] evaluate_0.14 BiocManager_1.30.16
#> [91] png_0.1-7 vctrs_0.3.8
#> [93] httpuv_1.6.3 tidyr_1.1.4
#> [95] gtable_0.3.0 purrr_0.3.4
#> [97] assertthat_0.2.1 cachem_1.0.6
#> [99] xfun_0.27 mime_0.12
#> [101] xtable_1.8-4 restfulr_0.0.13
#> [103] AnnotationFilter_1.18.0 later_1.3.0
#> [105] tibble_3.1.5 plyranges_1.14.0
#> [107] GenomicAlignments_1.30.0 AnnotationDbi_1.56.0
#> [109] memoise_2.0.0 IRanges_2.28.0
#> [111] ellipsis_0.3.2 interactiveDisplayBase_1.32.0