Contents

1 FLAMES

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 was created to simplify installation and execution of the FLAMES python module, which supports the article Comprehensive characterization of single cell full-length isoforms in human and mouse with long-read sequencing by Tian et al. (2020).

FLAMES pipeline workflow summary

Figure 1: FLAMES pipeline workflow summary

Input to FLAMES are fastq files generated from the long-read platform. Using a cell barcode allow-list (e.g. barcode list obtained from short-read sequencing), the find_barcode() function locates the barcode sequencing in the long-reads and trims the cell barcodes and the flanking UMI sequence. The pipeline then calls minimap_align(), a minimap2 wrapper, to align the demultiplex long-reads to the reference genome. Next, find_isoform() is called to identify novel isoforms using the alignment, creating an updated transcript assembly. The demultiplexed reads are then aligned to the transcript assembly by the function minimap_realign(). Finally, FLAMES calls quantify() 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.

2 FLAMES Pipeline Execution

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()).

2.1 FLAMES Single Cell Pipeline

2.1.1 Environment setup

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 match_barcode argument 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")
#> Warning: replacing previous import 'utils::findMatches' by
#> 'S4Vectors::findMatches' when loading 'AnnotationDbi'
#> Registered S3 method overwritten by 'GGally':
#>   method from   
#>   +.gg   ggplot2
#> Writing configuration parameters to:  /tmp/RtmpmQSwX4/filefaee7785c42ba/config_file_1027815.json

The 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_5end", min_sup_cnt = 10) with create a copy of the default config for 5’ 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.

2.1.2 FLAMES execution

Once the initial variables have been setup, the pipeline can be run using:

library(FLAMES)
# do not run if minimap2 cannot be found
if (is.character(locate_minimap2_dir())) { 
  sce <- sc_long_pipeline(
    annotation = annot, fastq = fastq, genome_fa = genome_fa,
    outdir = outdir, config_file = config_file, match_barcode = FALSE
  )
}
#> Warning in base::system2(command = "which", args = c("minimap2"), stdout =
#> TRUE, : running command ''which' minimap2 2>&1' had status 1

As stated above, the example dataset has already been demultiplexed, so match_barcode=FALSE should be set to skip the preprocessing of the barcodes.

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)
minimap2_dir <- locate_minimap2_dir()
#> Warning in base::system2(command = "which", args = c("minimap2"), stdout =
#> TRUE, : running command ''which' minimap2 2>&1' had status 1
if (is.character(minimap2_dir)) { # do not run if minimap2 cannot be found
  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, minimap2_dir = minimap2_dir
  ))
  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, minimap2_dir = minimap2_dir
  )
  quantify(annotation = annot, outdir = outdir, config = config)
  sce <- create_sce_from_dir(outdir = outdir, annotation = annot)
}

2.2 Visulisation

The combine_sce() and sc_annotate_plots() can be used for visulisation with experiment designs similar to FLT-seq. The combine_sce() function combines the SingleCellExperiment objects of 1) the short-read gene counts from the larger subsample, 2) the short-read gene counts from the smaller subsample and 3) the transcript counts from the larger subsample into a MultiAssayExperiment object. The sc_annotate_plots() function can then use it to produce annotated plots.

library(FLAMES)
library(MultiAssayExperiment)
 # load scm_lib80, scm_lib80 and scm_lib20
source(system.file("examples/scmixology1.R", package = "FLAMES"))
scm_lib80 <- scuttle::addPerCellQC(scm_lib80)
scm_lib20 <- scuttle::addPerCellQC(scm_lib20)
qc_80 <- scuttle::quickPerCellQC(colData(scm_lib80))
qc_20 <- scuttle::quickPerCellQC(colData(scm_lib20))
qc_80$discard[scm_lib80$sum < 10000] <- TRUE
qc_20$discard[scm_lib20$sum < 20000] <- TRUE

combined_sce <- combine_sce(
    short_read_large = scm_lib80[,!qc_80$discard],
    short_read_small = scm_lib20[,!qc_20$discard],
    long_read_sce = scm_lib20_transcripts,
    remove_duplicates = FALSE)

plots <- sc_annotate_plots(gene = "ENSG00000108107",
                           multiAssay = combined_sce)
#> Running PCA for experiments(multiAssay)$gene_counts ...
#> Running UMAP for experiments(multiAssay)$gene_counts ...
#> Imputing transcript counts ...
#> Plotting expression UMAPs ...
#> Cluster annotation not found, heatmaps skipped.

scater::plotReducedDim(
    experiments(plots[["multiAssay"]])$gene_counts,
    dimred = "PCA",
    colour_by = colData(plots[["multiAssay"]]
        [,,"gene_counts"])[,"Lib_small",drop=FALSE])


plots[["plot_umap"]]

plots[["combined_isoform_plot"]]

2.2.1 FLAMES termination

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 format
  • transcript_assembly.fa - transcript sequence from the isoforms
  • align2genome.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.

2.3 FLAMES Bulk Pipeline

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.

2.3.1 Environment setup

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:  /tmp/RtmpmQSwX4/filefaee7676fe405/config_file_1027815.json

2.3.2 FLAMES execution

Once the environment has been setup, the pipeline can be executed by running:

library(FLAMES)
if (is.character(locate_minimap2_dir())) {
  summarizedExperiment <- bulk_long_pipeline(
    annot = annot, fastq = fastq_dir, outdir = outdir,
    genome_fa = genome_fa, minimap2_dir = minimap2_dir,
    config_file = config_file
  )
}
#> Warning in base::system2(command = "which", args = c("minimap2"), stdout =
#> TRUE, : running command ''which' minimap2 2>&1' had status 1

2.3.3 FLAMES termination

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).

2.4 Running with Slurm

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_5end.json 
# set do_genome_alignment to true
sed -i '/do_genome_alignment/s/false/true/' \
        config_sclr_nanopore_5end.json 
srun -c 20 --mem 64GB \
        Rscript flames_pipeline.R && \

sed -i '/"do_/s/true/false/' \
        config_sclr_nanopore_5end.json && \
sed -i '/do_isoform_identification/s/false/true/' \
        config_sclr_nanopore_5end.json && \
srun -c 1 --mem 64GB --ntasks=1 \
        Rscript flames_pipeline.R  && \

sed -i '/"do_/s/true/false/' \
        config_sclr_nanopore_5end.json && \
sed -i '/do_read_realignment/s/false/true/' \
        config_sclr_nanopore_5end.json && \
srun -c 20 --mem 64GB --ntasks=1 \
        Rscript flames_pipeline.R && \

sed -i '/"do_/s/true/false/' \
        config_sclr_nanopore_5end.json && \
sed -i '/do_transcript_quantification/s/false/true/' \
        config_sclr_nanopore_5end.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.com

2.5 FLAMES on Windows

Due to FLAMES requiring minimap2 and pysam, FLAMES is currently unavaliable on Windows.

3 Session Info

#> 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] stats4    stats     graphics  grDevices utils     datasets  methods  
#> [8] base     
#> 
#> other attached packages:
#>  [1] SingleCellExperiment_1.22.0 MultiAssayExperiment_1.26.0
#>  [3] SummarizedExperiment_1.30.0 Biobase_2.60.0             
#>  [5] GenomicRanges_1.52.0        GenomeInfoDb_1.36.0        
#>  [7] IRanges_2.34.0              S4Vectors_0.38.0           
#>  [9] BiocGenerics_0.46.0         MatrixGenerics_1.12.0      
#> [11] matrixStats_0.63.0          FLAMES_1.6.0               
#> [13] BiocStyle_2.28.0           
#> 
#> loaded via a namespace (and not attached):
#>   [1] BiocIO_1.10.0             bitops_1.0-7             
#>   [3] filelock_1.0.2            tibble_3.2.1             
#>   [5] R.oo_1.25.0               basilisk.utils_1.12.0    
#>   [7] bambu_3.2.0               graph_1.78.0             
#>   [9] XML_3.99-0.14             rpart_4.1.19             
#>  [11] lifecycle_1.0.3           edgeR_3.42.0             
#>  [13] doParallel_1.0.17         OrganismDbi_1.42.0       
#>  [15] ensembldb_2.24.0          lattice_0.21-8           
#>  [17] backports_1.4.1           magrittr_2.0.3           
#>  [19] limma_3.56.0              Hmisc_5.0-1              
#>  [21] sass_0.4.5                rmarkdown_2.21           
#>  [23] jquerylib_0.1.4           yaml_2.3.7               
#>  [25] metapod_1.8.0             reticulate_1.28          
#>  [27] ggbio_1.48.0              cowplot_1.1.1            
#>  [29] DBI_1.1.3                 RColorBrewer_1.1-3       
#>  [31] zlibbioc_1.46.0           purrr_1.0.1              
#>  [33] R.utils_2.12.2            AnnotationFilter_1.24.0  
#>  [35] biovizBase_1.48.0         RCurl_1.98-1.12          
#>  [37] nnet_7.3-18               VariantAnnotation_1.46.0 
#>  [39] rappdirs_0.3.3            circlize_0.4.15          
#>  [41] GenomeInfoDbData_1.2.10   ggrepel_0.9.3            
#>  [43] irlba_2.3.5.1             dqrng_0.3.0              
#>  [45] DelayedMatrixStats_1.22.0 codetools_0.2-19         
#>  [47] DropletUtils_1.20.0       DelayedArray_0.26.0      
#>  [49] scuttle_1.10.0            xml2_1.3.3               
#>  [51] tidyselect_1.2.0          shape_1.4.6              
#>  [53] farver_2.1.1              viridis_0.6.2            
#>  [55] ScaledMatrix_1.8.0        BiocFileCache_2.8.0      
#>  [57] base64enc_0.1-3           GenomicAlignments_1.36.0 
#>  [59] jsonlite_1.8.4            BiocNeighbors_1.18.0     
#>  [61] GetoptLong_1.0.5          Formula_1.2-5            
#>  [63] scater_1.28.0             iterators_1.0.14         
#>  [65] foreach_1.5.2             tools_4.3.0              
#>  [67] progress_1.2.2            Rcpp_1.0.10              
#>  [69] glue_1.6.2                BiocBaseUtils_1.2.0      
#>  [71] gridExtra_2.3             xfun_0.39                
#>  [73] dplyr_1.1.2               HDF5Array_1.28.0         
#>  [75] withr_2.5.0               BiocManager_1.30.20      
#>  [77] fastmap_1.1.1             GGally_2.1.2             
#>  [79] basilisk_1.12.0           bluster_1.10.0           
#>  [81] rhdf5filters_1.12.0       fansi_1.0.4              
#>  [83] rsvd_1.0.5                digest_0.6.31            
#>  [85] R6_2.5.1                  colorspace_2.1-0         
#>  [87] dichromat_2.0-0.1         biomaRt_2.56.0           
#>  [89] RSQLite_2.3.1             R.methodsS3_1.8.2        
#>  [91] utf8_1.2.3                tidyr_1.3.0              
#>  [93] generics_0.1.3            data.table_1.14.8        
#>  [95] FNN_1.1.3.2               rtracklayer_1.60.0       
#>  [97] prettyunits_1.1.1         httr_1.4.5               
#>  [99] htmlwidgets_1.6.2         uwot_0.1.14              
#> [101] pkgconfig_2.0.3           gtable_0.3.3             
#> [103] blob_1.2.4                ComplexHeatmap_2.16.0    
#> [105] XVector_0.40.0            htmltools_0.5.5          
#> [107] bookdown_0.33             RBGL_1.76.0              
#> [109] ProtGenerics_1.32.0       clue_0.3-64              
#> [111] scales_1.2.1              png_0.1-8                
#> [113] scran_1.28.0              knitr_1.42               
#> [115] rstudioapi_0.14           tzdb_0.3.0               
#> [117] reshape2_1.4.4            rjson_0.2.21             
#> [119] checkmate_2.1.0           curl_5.0.0               
#> [121] cachem_1.0.7              rhdf5_2.44.0             
#> [123] GlobalOptions_0.1.2       stringr_1.5.0            
#> [125] vipor_0.4.5               parallel_4.3.0           
#> [127] foreign_0.8-84            AnnotationDbi_1.62.0     
#> [129] restfulr_0.0.15           GEOquery_2.68.0          
#> [131] pillar_1.9.0              grid_4.3.0               
#> [133] reshape_0.8.9             vctrs_0.6.2              
#> [135] BiocSingular_1.16.0       dbplyr_2.3.2             
#> [137] beachmat_2.16.0           cluster_2.1.4            
#> [139] beeswarm_0.4.0            htmlTable_2.4.1          
#> [141] evaluate_0.20             magick_2.7.4             
#> [143] readr_2.1.4               GenomicFeatures_1.52.0   
#> [145] cli_3.6.1                 locfit_1.5-9.7           
#> [147] compiler_4.3.0            Rsamtools_2.16.0         
#> [149] rlang_1.1.0               crayon_1.5.2             
#> [151] labeling_0.4.2            ggbeeswarm_0.7.1         
#> [153] plyr_1.8.8                stringi_1.7.12           
#> [155] viridisLite_0.4.1         BiocParallel_1.34.0      
#> [157] munsell_0.5.0             Biostrings_2.68.0        
#> [159] lazyeval_0.2.2            Matrix_1.5-4             
#> [161] dir.expiry_1.8.0          BSgenome_1.68.0          
#> [163] hms_1.1.3                 sparseMatrixStats_1.12.0 
#> [165] bit64_4.0.5               ggplot2_3.4.2            
#> [167] Rhdf5lib_1.22.0           statmod_1.5.0            
#> [169] KEGGREST_1.40.0           highr_0.10               
#> [171] igraph_1.4.2              memoise_2.0.1            
#> [173] bslib_0.4.2               bit_4.0.5                
#> [175] xgboost_1.7.5.1

References

Amezquita, 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.

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.