CellBarcode 1.6.0
library(data.table)
library(ggplot2)
library(CellBarcode)
The the information in the bam file. There are RNA sequence, Cell barcode and UMI in the bam file. We need to get the barcode in the RNA sequence together with the Cell barcode and UMI.
Where to find the bam file. Usually, the bam file in in following location of the CellRanger output.
CellRanger Output fold/outs/possorted_genome_bam.bam
Why preprocess. We need get the sam file as input. And in some cases we can do some filtering to make the input file smaller, for reducing the running time.
Example 1 get the sam file
samtools view possorted_genome_bam.bam > scRNASeq_10X.sam
Example 2 get the sam file only contain un mapped reads
In most of the time, the barcodes are designed not overlap with the genome sequence, so those barcodes sequences are not mapped to the reference genome. Add a simple parameter, we can get the un-mapped reads to significantly reduce the running time of the barcode extraction procedure. In the following example, the scRNASeq_10X.sam
file only contains the un-mapped reads.
samtools view -f 4 possorted_genome_bam.bam > scRNASeq_10X.sam
Extract lineage barcode with cell-barcode, UMI information.
The parameters:
sam: The location of sam format file derived from CellRanger output bam file.
pattern: A regular expression, that matchs the lineage barcode. The first capture will be extract as lineage barcode. Please check next section to know more about how to define a pattern.
cell_barcode_tag: The cell barcode field tag in the sam file, the default is “CR” in the Cell Ranger output. You do not need to modify it, unless you know what you are doing.
umi_tag: The UMI field tag in the sam file, the default is “UR” in the Cell Ranger output. You do not need to modify it, unless you know what you are doing.
sam_file <- system.file("extdata", "scRNASeq_10X.sam", package = "CellBarcode")
d = bc_extract_10XscSeq(
sam = sam_file,
pattern = "AGATCAG(.*)TGTGGTA",
cell_barcode_tag = "CR",
umi_tag = "UR"
)
head(d)
#> cell_barcode umi barcode_seq count
#> 1 AATAGAGGTTCTAACG AAGACAGGCTTA AGCCATTATGTCGTGTCATGA 1
#> 2 AAGTGAAAGTCCCGGT AGCCATTGTCCT TTTAATAGTTTCTATTTTGGA 1
#> 3 AAGTCGTAGCACGGAT AAAGCCTGCAAG CTATATATTCGTTGTATTATA 1
#> 4 AAGACAAAGGCGCTTC CGTCTTGGTAGC TCATTTGATACTAAGAATGAA 2
#> 5 AAGACAAAGCTAAATG GTATTCTCTGAG ATGGTTGTAGATG 1
#> 6 AACTTCTGTGTTCCTC GACAACCTTGGT GGGTGGGATTCAGTGATCATA 1
The output is a data.frame
. It contains 3 columns:
cell_barcode: It is unique identificaiton for each cell. One cell_barcode
corresponding to a cell.
umi: It is used to label the molecular before PCR, a unique UMI means one molecular.
barcode_seq: Depends on the experiment design, it usually labels the lineage.
count: The reads number of a specifc cell-barcode, UMI, lineage barcode combination.
pattern
The pattern
is a regular expression, it tells the function where
to find the barcode. In the pattern, we define the barcode backbone, and label the barcode sequence by bracket ()
.
For example, the pattern ATCG(.{21})TCGG
tells the barcode is surrounded by constant sequence of ATCG
, and TCGG
. Following are some examples to define the constant region and barcode sequence.
Example 1
ATCG(.{21})
21 bases barcode after a constant sequence of “ATCG”.
Example 2
(.{15})TCGA
15 bases barcode before a constant sequence of “TCGA”.
Example 3
ATCG(.*)TCGA
A flexible length barcode between constant regions of “ATCG” and “TCGA”.
Need more helps: About more complex barcode pattern, please ask the package author, then the exmaple will be apear here.
sessionInfo()
#> R version 4.3.0 RC (2023-04-13 r84269)
#> Platform: x86_64-pc-linux-gnu (64-bit)
#> Running under: Ubuntu 22.04.2 LTS
#>
#> Matrix products: default
#> BLAS: /home/biocbuild/bbs-3.17-bioc/R/lib/libRblas.so
#> LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.10.0
#>
#> 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
#>
#> time zone: America/New_York
#> tzcode source: system (glibc)
#>
#> attached base packages:
#> [1] stats graphics grDevices utils datasets methods base
#>
#> other attached packages:
#> [1] CellBarcode_1.6.0 ggplot2_3.4.2 data.table_1.14.8 BiocStyle_2.28.0
#>
#> loaded via a namespace (and not attached):
#> [1] SummarizedExperiment_1.30.0 gtable_0.3.3
#> [3] xfun_0.39 bslib_0.4.2
#> [5] hwriter_1.3.2.1 latticeExtra_0.6-30
#> [7] Biobase_2.60.0 lattice_0.21-8
#> [9] Rdpack_2.4 vctrs_0.6.2
#> [11] tools_4.3.0 bitops_1.0-7
#> [13] generics_0.1.3 stats4_4.3.0
#> [15] parallel_4.3.0 tibble_3.2.1
#> [17] fansi_1.0.4 pkgconfig_2.0.3
#> [19] Matrix_1.5-4 RColorBrewer_1.1-3
#> [21] S4Vectors_0.38.0 lifecycle_1.0.3
#> [23] GenomeInfoDbData_1.2.10 stringr_1.5.0
#> [25] deldir_1.0-6 compiler_4.3.0
#> [27] egg_0.4.5 Rsamtools_2.16.0
#> [29] Biostrings_2.68.0 Ckmeans.1d.dp_4.3.4
#> [31] munsell_0.5.0 codetools_0.2-19
#> [33] GenomeInfoDb_1.36.0 htmltools_0.5.5
#> [35] sass_0.4.5 RCurl_1.98-1.12
#> [37] yaml_2.3.7 pillar_1.9.0
#> [39] crayon_1.5.2 jquerylib_0.1.4
#> [41] BiocParallel_1.34.0 cachem_1.0.7
#> [43] DelayedArray_0.26.0 ShortRead_1.58.0
#> [45] tidyselect_1.2.0 digest_0.6.31
#> [47] stringi_1.7.12 dplyr_1.1.2
#> [49] bookdown_0.33 fastmap_1.1.1
#> [51] grid_4.3.0 colorspace_2.1-0
#> [53] cli_3.6.1 magrittr_2.0.3
#> [55] utf8_1.2.3 withr_2.5.0
#> [57] scales_1.2.1 rmarkdown_2.21
#> [59] XVector_0.40.0 jpeg_0.1-10
#> [61] matrixStats_0.63.0 interp_1.1-4
#> [63] gridExtra_2.3 png_0.1-8
#> [65] evaluate_0.20 knitr_1.42
#> [67] rbibutils_2.2.13 GenomicRanges_1.52.0
#> [69] IRanges_2.34.0 rlang_1.1.0
#> [71] Rcpp_1.0.10 glue_1.6.2
#> [73] BiocManager_1.30.20 BiocGenerics_0.46.0
#> [75] jsonlite_1.8.4 plyr_1.8.8
#> [77] R6_2.5.1 MatrixGenerics_1.12.0
#> [79] GenomicAlignments_1.36.0 zlibbioc_1.46.0