FLAMES 2.0.2
The FLAMES package provides a framework for performing single-cell and bulk read full-length analysis of mutations and splicing. FLAMES performs cell barcode and UMI assignment from nanopore reads as well as semi-supervised isoform detection and quantification. FLAMES is designed to be an easy and quick to use, powerful workflow for isoform detection and quantification, splicing analysis and mutation detection, and seeks to overcome the limitations of other tools, such as an inability to process single cell data, and a focus on cell barcode and UMI assignment (Tian et al. 2020).
This R package represents an enhanced iteration of the FLAMES Python module that originally designed to support the research work presented in Comprehensive characterization of single-cell full-length isoforms in human and mouse with long-read sequencing by Tian et al (2020). This upgraded R version not only simplifies the installation and execution processes but also incorporates additional functionality, streamlining the analysis of single-cell full-length isoforms using long-read sequencing data.
Figure 1: FLAMES pipeline
When processing single cell data, FLAMES can be run on data generated from long-read platform with or without matched short-read sequencing data. If only the long reads are available, FLAMES takes as input the long reads in fastq files. The blaze() function is then called to demultiplex the long reads into cell-specific fastq files. If the short-read data are unavailable, FLAMES incorporates the BLAZE as blaze() function to identify cell barcodes solely from long reads (You et al. 2023). The blaze() function is called to locate the barcode sequencing in the long-reads, identify barcode allow-list directly from long reads, assign reads to cells (i.e., demultiplexing) and trims the cell barcodes and the flanking UMI sequence.
When matched short-read single-cell RNA sequencing data is available, FLAMES could take the cell barcode allow-list generated from short reads in addition to the long-read fastq files. The short-read allow-list will used as a reference to guide the demultiplexing of the long reads. To do this, FLAMES incorporates the flexiplex as find_barcode() function to perform demultiplexing (Davidson et al. 2023). The find_barcode() function is called to extract the cell barcode and UMI sequence from the long reads, using the allow-list as a reference.
After the demultiplexing, the pipeline calls minimap_align(), a minimap2 wrapper, to align the demultiplex long-reads to the reference genome. quantify_gene() is then called to deduplicate reads from same unique molecular identifiers (UMIs) and generate gene-level UMI counts. Note that during the UMI deduplication, the pipeline only keeps the longest reads among those with the same UMI for downstream steps.
Next, find_isoform() is called to identify novel isoforms using the alignment, creating an updated transcript assembly. In this step, the users may choose to use bambu, which is designed for transcript discovery and quantification using long read RNA-Seq data (Chen et al. 2023). Otherwise, users could use the build-in isoform identification methods in FLAMES. Both methods have been wrapped in the find_isoform() function.
Afterwards, the demultiplexed reads are aligned to the transcript assembly by the function minimap_realign(). Finally, FLAMES calls quantify_transcript() to quantify the transcript counts using the (re-)alignment file.
Figure 1 summerise the steps of the pipeline described above. The pipeline functions (sc_long_pipeline, bulk_long_pipeline, sc_long_multisample_pipline) execute the steps sequentially and return SingleCellExperiment object, SummarizedExperiment object and a list of SingleCellExperiment objects respectivly.
For read alignment and realignment, FLAMES uses minimap2, a versatile alignment program for aligning sequences against a reference database (Li 2018). This software needs to be downloaded prior to using the FLAMES pipeline, and can be found at https://github.com/lh3/minimap2.
This vignette will detail the process of running the FLAMES pipeline. It details the execution of both the single cell pipeline (sc_long_pipeline()) and the bulk data pipeline (bulk_long_pipeline()).
To get started, the pipeline needs access to a gene annotation file in GFF3 or GTF format, a directory containing one or more FASTQ files (which will be merged as pre-processing), a genome FASTA file, as well as the file path to minimap2, and the file path to the directory to hold output files.
The single cell pipeline can demultiplex the input data before running, if required. This can be disabled by setting the do_barcode_demultiplex argument in the config file to FALSE when calling the pipeline. This example dataset has already been demultiplexed.
For this example, the required files are downloaded from GitHub using BiocFileCache (Shepherd and Morgan 2021).
temp_path <- tempfile()
bfc <- BiocFileCache::BiocFileCache(temp_path, ask = FALSE)
file_url <-
  "https://raw.githubusercontent.com/OliverVoogd/FLAMESData/master/data"
annot <- bfc[[names(BiocFileCache::bfcadd(
  bfc, "Annotation",
  file.path(file_url, "gencodeshortened.gtf")
))]]
genome_fa <- bfc[[names(BiocFileCache::bfcadd(
  bfc,
  "Genomefa",
  file.path(file_url, "GRCh38shortened.fa")
))]]
fastq <- bfc[[names(BiocFileCache::bfcadd(
  bfc, "Fastq", file.path(file_url, "sc_align2genome.sample.fastq.gz")))]]
# setup other environment variables
outdir <- tempfile()
dir.create(outdir)
config_file <- FLAMES::create_config(outdir, type = "SIRV", do_barcode_demultiplex = TRUE)
#> Registered S3 method overwritten by 'GGally':
#>   method from   
#>   +.gg   ggplot2
#> Writing configuration parameters to:  /home/biocbuild/bbs-3.20-bioc/tmpdir/Rtmp2fZaYC/file1bd7639aac5a/config_file_114038.jsonThe optional argument config_file can be given to both bulk_long_pipeline() and sc_long_pipeline() in order to customise the execution of the pipelines. It is expected to be a JSON file, and an example can be found by calling FLAMES::create_config, which returns the path to a copy of the default JSON configuration file. The values from the default configs can be altered by either editing the JSON file manually, or passing additional arguments to FLAMES::create_config, for example, FLAMES::create_config(outdir, type = "sc_3end", min_sup_cnt = 10) with create a copy of the default config for 10X 3’ end nanopore single cell reads, with the min_sup_cnt value (minimal number of supporting reads for a novel isoform to pass filtering) changed to 10.
This vignette uses the default configuration file created for SIRV reads.
Once the initial variables have been setup, the pipeline can be run using:
library(FLAMES)
# do not run if minimap2 cannot be found
if (!any(is.na(find_bin(c("minimap2", "k8"))))) {
  sce <- sc_long_pipeline(
    annotation = annot, fastq = fastq, genome_fa = genome_fa,
    outdir = outdir, config_file = config_file, expect_cell_number = 10)
}If, however, the input fastq files need to be demultiplexed, then the reference_csv argument will need to be specified - reference_csv is the file path to a cell barcode allow-list, as a text file with one barcode per line. The filtered_feature_bc_matrix/barcodes.tsv.gz from cellranger outputs can be used to create such allow-list, with zcat barcodes.tsv.gz | cut -f1 -d'-' > allow-list.csv.
The pipeline can alse be run by calling the consituent steps sequentially:
library(FLAMES)
# do not run if minimap2 cannot be found
if (!any(is.na(find_bin(c("minimap2", "k8"))))) {
  config <- jsonlite::fromJSON(config_file)
  # find_barcode(...)
  genome_bam <- rownames(minimap2_align(
    config = config, fa_file = genome_fa, fq_in = fastq, annot = annot,
    outdir = outdir
  ))
  find_isoform(
    annotation = annot, genome_fa = genome_fa,
    genome_bam = genome_bam, outdir = outdir, config = config
  )
  minimap2_realign(
    config = config, fq_in = fastq,
    outdir = outdir
  )
  quantify_transcript(annotation = annot, outdir = outdir, config = config)
  sce <- create_sce_from_dir(outdir = outdir, annotation = annot)
}The plot_isoform_reduced_dim() function can be used to visualise the isoform expression in reduced dimensions. The function takes a SingleCellExperiment object and either a gene ID or a vector of transcript IDs as input. If the SingleCellExperiment object contains gene counts with the transcript counts stored as the “transcript” altExp slot, the function will use the reduced dimensions of the gene counts and plot the transcript expression values.
library(FLAMES)
library(SingleCellExperiment)
# just the transcript counts
scmixology_lib10_transcripts |>
  scuttle::logNormCounts() |>
  scater::runPCA() |>
  scater::runUMAP() |>
  plot_isoform_reduced_dim('ENSG00000108107')
# SCE with gene counts as main assay and transcript counts as altExp
scmixology_lib10 <- scmixology_lib10[, colSums(counts(scmixology_lib10)) > 0]
sce_lr <- scmixology_lib10[, colnames(scmixology_lib10) %in% colnames(scmixology_lib10_transcripts)]
altExp(sce_lr, "transcript") <- scmixology_lib10_transcripts[, colnames(sce_lr)]
sce_lr |>
  scuttle::logNormCounts() |>
  scater::runPCA() |>
  scater::runUMAP() |>
  plot_isoform_reduced_dim('ENSG00000108107')For experiments where the long-read is preformed on a subsample of cells and short-reads are preformed on the full sample, the combine_sce() function can be used to combine the SingleCellExperiment objects of the two samples. The first argument should be the subsample where the main assay is the gene counts and the transcript counts is saved as the “transcript” altExp slot. The second argument should be the other sample with only gene counts. This function returns a SingleCellExperiment object with the main assay containing the combined gene counts and a “transcript” altExp slot, where cells missing long-read data have NA values.
The sc_impute_transcript() function can be then used to impute the missing transcript counts using the gene counts.
The plot_isoform_reduced_dim() function can handle the SingleCellExperiment object with the combined gene and transcript counts, with or without imputed values, NA values will be plotted in grey.
combined_sce <- combine_sce(sce_lr, scmixology_lib90)
combined_sce <- combined_sce |>
  scuttle::logNormCounts() |>
  scater::runPCA() |>
  scater::runUMAP()
# plot without imputation
plot_isoform_reduced_dim(combined_sce, 'ENSG00000108107')
# impute missing transcript counts
combined_sce_impute <- sc_impute_transcript(combined_sce)
#> Imputing transcript counts ...
plot_isoform_reduced_dim(combined_sce_impute, 'ENSG00000108107')The directory outdir now contains several output files returned from this pipeline. The output files generated by this pipeline are:
transcript_count.csv.gz - a transcript count matrix (also contained in the output SummarizedExperiment or SingleCellExperiment)isoform_annotated.filtered.gff3 - found isoform information in gff3 formattranscript_assembly.fa - transcript sequence from the isoformsalign2genome.bam sorted BAM file with reads aligned to genome (intermediate FLAMES step)realign2transcript.bam - sorted realigned BAM file using the transcript_assembly.fa as reference (intermediate FLAMES step)tss_tes.bedgraph- TSS TES enrichment for all reads (for QC)The pipeline also returns a SummarizedExperiment or SingleCellExperiment object, depending on the pipeline run, containing the data from transcript_count.csv.gzand isoform_annotated.filtered.gff3 (Amezquita et al. 2020) (Morgan et al. 2020). This SummarizedExperiment (or SingleCellExperiment) object contains the same data as present in the outdir directory, and is given to simplify the process of reading the transcript count matrix and annotation data back into R, for further analysis.
A basic example of the execution of the FLAMES bulk pipeline is given below. The process for this is essentially identical to the above example for single cell data.
To get started, the pipeline needs access to a gene annotation file in GFF3 or GTF format, a directory containing one or more FASTQ files (which will be merged as pre-processing), a genome FASTA file, as well as the file path to minimap2, and the file path to the directory to hold output files.
For this example, these files are downloaded from GitHub using BiocFileCache (Shepherd and Morgan 2021).
temp_path <- tempfile()
bfc <- BiocFileCache::BiocFileCache(temp_path, ask = FALSE)
file_url <-
  "https://raw.githubusercontent.com/OliverVoogd/FLAMESData/master/data"
annot <- bfc[[names(BiocFileCache::bfcadd(
  bfc, "Annotation",
  file.path(file_url, "SIRV_isoforms_multi-fasta-annotation_C_170612a.gtf")
))]]
genome_fa <- bfc[[names(BiocFileCache::bfcadd(
  bfc, "Genomefa",
  file.path(file_url, "SIRV_isoforms_multi-fasta_170612a.fasta")
))]]
# download the two fastq files, move them to a folder to be merged together
fastq1 <- bfc[[names(BiocFileCache::bfcadd(bfc, "Fastq1", 
                file.path(file_url, "fastq", "sample1.fastq.gz")))]]
fastq2 <- bfc[[names(BiocFileCache::bfcadd(bfc, "Fastq2",
                file.path(file_url, "fastq", "sample2.fastq.gz")))]]
# the downloaded fastq files need to be in a directory to be merged together
fastq_dir <- file.path(temp_path, "fastq_dir") 
dir.create(fastq_dir)
file.copy(c(fastq1, fastq2), fastq_dir)
#> [1] TRUE TRUE
unlink(c(fastq1, fastq2)) # the original files can be deleted
# setup other environment variables
outdir <- tempfile()
dir.create(outdir)
config_file <- FLAMES::create_config(outdir)
#> Writing configuration parameters to:  /home/biocbuild/bbs-3.20-bioc/tmpdir/Rtmp2fZaYC/file1bd762e9d807b/config_file_114038.jsonOnce the environment has been setup, the pipeline can be executed by running:
library(FLAMES)
if (!any(is.na(find_bin(c("minimap2", "k8"))))) {
  summarizedExperiment <- bulk_long_pipeline(
    annot = annot, fastq = fastq_dir, outdir = outdir,
    genome_fa = genome_fa, config_file = config_file
  )
}After the bulk pipeline has completed, the output directory contains the same files as the single cell pipeline produces. bulk_long_pipeline also returns a SummarizedExperiment object, containing the same data as the SingleCellExperiment as above (Amezquita et al. 2020) (Morgan et al. 2020).
Since the barcode demultiplexing step, isoform identification step and quantification step is currently single threaded, job scheduler (such as Slurm) users may want to allocate different resources for each step. This can be achieved by calling the individual functions (find_isoform, minimap2_realign etc.) sequentially, or by changing the configuration file. The Bash script below demonstrates how this can be done with a single R script calling the pipeline function.
#!/bin/bash -x
# set all steps to false
sed -i '/"do_/s/true/false/' config_sclr_nanopore_3end.json 
# set do_genome_alignment to true
sed -i '/do_genome_alignment/s/false/true/' \
        config_sclr_nanopore_3end.json 
srun -c 20 --mem 64GB \
        Rscript flames_pipeline.R && \
sed -i '/"do_/s/true/false/' \
        config_sclr_nanopore_3end.json && \
sed -i '/do_isoform_identification/s/false/true/' \
        config_sclr_nanopore_3end.json && \
srun -c 1 --mem 64GB --ntasks=1 \
        Rscript flames_pipeline.R  && \
sed -i '/"do_/s/true/false/' \
        config_sclr_nanopore_3end.json && \
sed -i '/do_read_realignment/s/false/true/' \
        config_sclr_nanopore_3end.json && \
srun -c 20 --mem 64GB --ntasks=1 \
        Rscript flames_pipeline.R && \
sed -i '/"do_/s/true/false/' \
        config_sclr_nanopore_3end.json && \
sed -i '/do_transcript_quantification/s/false/true/' \
        config_sclr_nanopore_3end.json && \
srun -c 1 --mem 64GB --ntasks=1 \
        Rscript flames_pipeline.R && \
echo "Pipeline finished"The Bash script can then be executed inside a tmux or screen session with:
./flames.sh | tee -a bash.log; tail bash.log \
        | mail -s "flames pipeline ended" user@example.comDue to FLAMES requiring minimap2 and pysam, FLAMES is currently unavaliable on Windows.
Please cite flames’s paper if you use flames in your research. As FLAMES used incorporates BLAZE, flexiplex and minimap2, samtools, bambu. Please make sure to cite when using these tools.
#> R version 4.4.2 (2024-10-31)
#> Platform: x86_64-pc-linux-gnu
#> Running under: Ubuntu 24.04.1 LTS
#> 
#> Matrix products: default
#> BLAS:   /home/biocbuild/bbs-3.20-bioc/R/lib/libRblas.so 
#> LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.12.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] stats4    stats     graphics  grDevices utils     datasets  methods  
#> [8] base     
#> 
#> other attached packages:
#>  [1] SingleCellExperiment_1.28.1 SummarizedExperiment_1.36.0
#>  [3] Biobase_2.66.0              GenomicRanges_1.58.0       
#>  [5] GenomeInfoDb_1.42.1         IRanges_2.40.1             
#>  [7] S4Vectors_0.44.0            BiocGenerics_0.52.0        
#>  [9] MatrixGenerics_1.18.0       matrixStats_1.4.1          
#> [11] FLAMES_2.0.2                BiocStyle_2.34.0           
#> 
#> loaded via a namespace (and not attached):
#>   [1] fs_1.6.5                  ProtGenerics_1.38.0      
#>   [3] bitops_1.0-9              httr_1.4.7               
#>   [5] RColorBrewer_1.1-3        doParallel_1.0.17        
#>   [7] tools_4.4.2               backports_1.5.0          
#>   [9] utf8_1.2.4                R6_2.5.1                 
#>  [11] HDF5Array_1.34.0          uwot_0.2.2               
#>  [13] lazyeval_0.2.2            rhdf5filters_1.18.0      
#>  [15] GetoptLong_1.0.5          withr_3.0.2              
#>  [17] prettyunits_1.2.0         GGally_2.2.1             
#>  [19] gridExtra_2.3             cli_3.6.3                
#>  [21] scatterpie_0.2.4          labeling_0.4.3           
#>  [23] ggbio_1.54.0              sass_0.4.9               
#>  [25] readr_2.1.5               Rsamtools_2.22.0         
#>  [27] yulab.utils_0.1.8         txdbmaker_1.2.1          
#>  [29] foreign_0.8-87            R.utils_2.12.3           
#>  [31] dichromat_2.0-0.1         scater_1.34.0            
#>  [33] parallelly_1.40.1         BSgenome_1.74.0          
#>  [35] limma_3.62.1              rstudioapi_0.17.1        
#>  [37] RSQLite_2.3.9             FNN_1.1.4.1              
#>  [39] generics_0.1.3            shape_1.4.6.1            
#>  [41] BiocIO_1.16.0             dplyr_1.1.4              
#>  [43] Matrix_1.7-1              ggbeeswarm_0.7.2         
#>  [45] fansi_1.0.6               abind_1.4-8              
#>  [47] R.methodsS3_1.8.2         lifecycle_1.0.4          
#>  [49] yaml_2.3.10               edgeR_4.4.1              
#>  [51] rhdf5_2.50.1              SparseArray_1.6.0        
#>  [53] BiocFileCache_2.14.0      grid_4.4.2               
#>  [55] blob_1.2.4                dqrng_0.4.1              
#>  [57] crayon_1.5.3              dir.expiry_1.14.0        
#>  [59] lattice_0.22-6            beachmat_2.22.0          
#>  [61] cowplot_1.1.3             GenomicFeatures_1.58.0   
#>  [63] KEGGREST_1.46.0           magick_2.8.5             
#>  [65] metapod_1.14.0            pillar_1.9.0             
#>  [67] knitr_1.49                ComplexHeatmap_2.22.0    
#>  [69] bambu_3.8.0               rjson_0.2.23             
#>  [71] xgboost_1.7.8.1           codetools_0.2-20         
#>  [73] glue_1.8.0                ggfun_0.1.8              
#>  [75] data.table_1.16.4         vctrs_0.6.5              
#>  [77] png_0.1-8                 gtable_0.3.6             
#>  [79] cachem_1.1.0              xfun_0.49                
#>  [81] S4Arrays_1.6.0            DropletUtils_1.26.0      
#>  [83] iterators_1.0.14          tinytex_0.54             
#>  [85] bluster_1.16.0            statmod_1.5.0            
#>  [87] bit64_4.5.2               progress_1.2.3           
#>  [89] filelock_1.0.3            bslib_0.8.0              
#>  [91] irlba_2.3.5.1             vipor_0.4.7              
#>  [93] rpart_4.1.23              colorspace_2.1-1         
#>  [95] DBI_1.2.3                 Hmisc_5.2-1              
#>  [97] nnet_7.3-19               tidyselect_1.2.1         
#>  [99] bit_4.5.0.1               compiler_4.4.2           
#> [101] curl_6.0.1                httr2_1.0.7              
#> [103] graph_1.84.0              htmlTable_2.4.3          
#> [105] BiocNeighbors_2.0.1       basilisk.utils_1.18.0    
#> [107] xml2_1.3.6                DelayedArray_0.32.0      
#> [109] bookdown_0.41             rtracklayer_1.66.0       
#> [111] checkmate_2.3.2           scales_1.3.0             
#> [113] RBGL_1.82.0               rappdirs_0.3.3           
#> [115] stringr_1.5.1             SpatialExperiment_1.16.0 
#> [117] digest_0.6.37             rmarkdown_2.29           
#> [119] basilisk_1.18.0           XVector_0.46.0           
#> [121] htmltools_0.5.8.1         pkgconfig_2.0.3          
#> [123] base64enc_0.1-3           sparseMatrixStats_1.18.0 
#> [125] dbplyr_2.5.0              fastmap_1.2.0            
#> [127] ensembldb_2.30.0          rlang_1.1.4              
#> [129] GlobalOptions_0.1.2       htmlwidgets_1.6.4        
#> [131] UCSC.utils_1.2.0          DelayedMatrixStats_1.28.0
#> [133] farver_2.1.2              jquerylib_0.1.4          
#> [135] jsonlite_1.8.9            BiocParallel_1.40.0      
#> [137] R.oo_1.27.0               BiocSingular_1.22.0      
#> [139] VariantAnnotation_1.52.0  RCurl_1.98-1.16          
#> [141] magrittr_2.0.3            Formula_1.2-5            
#> [143] scuttle_1.16.0            GenomeInfoDbData_1.2.13  
#> [145] Rhdf5lib_1.28.0           munsell_0.5.1            
#> [147] Rcpp_1.0.13-1             viridis_0.6.5            
#> [149] reticulate_1.40.0         stringi_1.8.4            
#> [151] zlibbioc_1.52.0           MASS_7.3-61              
#> [153] plyr_1.8.9                ggstats_0.7.0            
#> [155] parallel_4.4.2            listenv_0.9.1            
#> [157] ggrepel_0.9.6             Biostrings_2.74.0        
#> [159] hms_1.1.3                 circlize_0.4.16          
#> [161] locfit_1.5-9.10           igraph_2.1.2             
#> [163] reshape2_1.4.4            biomaRt_2.62.0           
#> [165] ScaledMatrix_1.14.0       XML_3.99-0.17            
#> [167] evaluate_1.0.1            biovizBase_1.54.0        
#> [169] scran_1.34.0              BiocManager_1.30.25      
#> [171] tzdb_0.4.0                foreach_1.5.2            
#> [173] tweenr_2.0.3              tidyr_1.3.1              
#> [175] purrr_1.0.2               polyclip_1.10-7          
#> [177] future_1.34.0             clue_0.3-66              
#> [179] ggplot2_3.5.1             ggforce_0.4.2            
#> [181] rsvd_1.0.5                restfulr_0.0.15          
#> [183] AnnotationFilter_1.30.0   viridisLite_0.4.2        
#> [185] OrganismDbi_1.48.0        tibble_3.2.1             
#> [187] memoise_2.0.1             beeswarm_0.4.0           
#> [189] AnnotationDbi_1.68.0      GenomicAlignments_1.42.0 
#> [191] cluster_2.1.7             globals_0.16.3Amezquita, Robert, Aaron Lun, Etienne Becht, Vince Carey, Lindsay Carpp, Ludwig Geistlinger, Federico Marini, et al. 2020. “Orchestrating Single-Cell Analysis with Bioconductor.” Nature Methods 17: 137–45. https://www.nature.com/articles/s41592-019-0654-x.
Chen, Ying, Andre Sim, Yuk Kei Wan, Keith Yeo, Joseph Jing Xian Lee, Min Hao Ling, Michael I Love, and Jonathan Göke. 2023. “Context-Aware Transcript Quantification from Long-Read Rna-Seq Data with Bambu.” Nature Methods, 1–9.
Davidson, Nadia M, Noorul Amin, Ling Min Hao, Changqing Wang, Oliver Cheng, Jonathan M Goeke, Matthew E Ritchie, and Shuyi Wu. 2023. “Flexiplex: A Versatile Demultiplexer and Search Tool for Omics Data.” bioRxiv, 2023–08.
Li, Heng. 2018. “Minimap2: pairwise alignment for nucleotide sequences.” Bioinformatics 34 (18): 3094–3100. https://doi.org/10.1093/bioinformatics/bty191.
Morgan, Martin, Valerie Obenchain, Jim Hester, and Hervé Pagès. 2020. SummarizedExperiment: SummarizedExperiment Container. https://bioconductor.org/packages/SummarizedExperiment.
Shepherd, Lori, and Martin Morgan. 2021. BiocFileCache: Manage Files Across Sessions.
Tian, Luyi, Jafar S. Jabbari, Rachel Thijssen, Quentin Gouil, Shanika L. Amarasinghe, Hasaru Kariyawasam, Shian Su, et al. 2020. “Comprehensive Characterization of Single Cell Full-Length Isoforms in Human and Mouse with Long-Read Sequencing.” bioRxiv. https://doi.org/10.1101/2020.08.10.243543.
You, Yupei, Yair DJ Prawer, De Paoli-Iseppi, Cameron PJ Hunt, Clare L Parish, Heejung Shim, Michael B Clark, and others. 2023. “Identification of Cell Barcodes from Long-Read Single-Cell Rna-Seq with Blaze.” Genome Biology 24 (1): 1–23.