--- title: "Analysing Long Read RNA-Seq data with bambu" author: Ying Chen, Andre Sim, Yuk Kei Wan, Jonathan Göke output: rmarkdown::html_vignette vignette: > %\VignetteIndexEntry{bambu} %\VignetteEngine{knitr::rmarkdown} %\VignettePackage{bambu} %\VignetteEncoding{UTF-8} --- ```{r, include = FALSE} knitr::opts_chunk$set( collapse = TRUE,tidy = TRUE, warning=FALSE, message=FALSE, comment = "##>" ) ``` # Introduction *[Bambu](https://github.com/GoekeLab/bambu)* is a method for transcript discovery and quantification from long read RNA-Seq data. Bambu uses aligned reads and genome reference annotations as input, and will return abundance estimates for all known transcripts and for newly discovered transcripts. Bambu uses the information from the reference annotations to correct misalignment at splice junctions, then reduces the aligned reads to read equivalent classes, and uses this information to identify novel transcripts across all samples of interest. Reads are then assigned to transcripts, and expression estimates are obtained using an expectation maximisation algorithm. Here, we present an example workflow for analysing Nanopore long read RNA-Sequencing data from two human cancer cell lines from the Singapore Nanopore Expression Project (SG-NEx). # Quick start: Transcript discovery and quantification with bambu ## Installation You can install bambu from github: ```{r, eval = FALSE} if (!requireNamespace("BiocManager", quietly = TRUE)) install.packages("BiocManager") BiocManager::install("bambu") BiocManager::install("NanoporeRNASeq") ``` ## General Usage The default mode to run *bambu* is using a set of aligned reads (bam files), reference genome annotations (gtf file, TxDb object, or bambuAnnotation object), and reference genome sequence (fasta file or BSgenome). *bambu* will return a summarizedExperiment object with the genomic coordinates for annotated and new transcripts and transcript expression estimates. We highly recommend to use the same annotations that were used for genome alignment. If you have a gtf file and fasta file you can run bambu with the following options: ```{r, results = "hide"} library(bambu) test.bam <- system.file("extdata", "SGNex_A549_directRNA_replicate5_run1_chr9_1_1000000.bam", package = "bambu") fa.file <- system.file("extdata", "Homo_sapiens.GRCh38.dna_sm.primary_assembly_chr9_1_1000000.fa", package = "bambu") gtf.file <- system.file("extdata", "Homo_sapiens.GRCh38.91_chr9_1_1000000.gtf", package = "bambu") bambuAnnotations <- prepareAnnotations(gtf.file) se <- bambu(reads = test.bam, annotations = bambuAnnotations, genome = fa.file) ``` *bambu* returns a SummarizedExperiment object which can be accessed as follows: * assays(se) returns the transcript abundance estimates as counts or CPM, as well the full length or unique read counts * rowRanges(se) returns a GRangesList with all annotated and newly discovered transcripts * rowData(se) returns additional information about each transcript such as the gene name and the class of newly discovered transcript # A complete workflow to identify and quantify transcript expression from Nanopore RNA-Seq data To demonstrate the usage of Bambu, we used long-read RNA-Seq data generated using Oxford Nanopore Sequencing from the NanoporeRNASeq package, which consists of 6 samples from two human cell lines (K562 and MCF7) that were generated by the [SG-NEx project](https://github.com/GoekeLab/sg-nex-data). Each of these cell lines has three replicates, with 1 direct RNA sequencing run and 2 cDNA sequencing runs. Reads are aligned to chromosome 22 (Grch38) and stored as bam files. In this workflow, we will demonstrate how to apply *bambu* to these bam files to identify novel transcripts and estimate transcript expression levels, visualize the results, and identify differentially expressed genes and transcripts. ## Input data ### Aligned reads (bam files) *bambu* takes genomic alignments saved in bam files. Here we use bam-files from the *NanoporeRNASeq* package, which contains reads aligned to the first half of the human chromosome 22 using *minimap2*. ```{r} library(bambu) library(NanoporeRNASeq) data("SGNexSamples") SGNexSamples library(ExperimentHub) NanoporeData <- query(ExperimentHub(), c("NanoporeRNA", "GRCh38","BAM")) bamFile <- Rsamtools::BamFileList(NanoporeData[["EH3808"]]) ``` ### Genome sequence (fasta file/ BSGenome object) *bambu* additionally requires a genome sequence, which is used to correct splicing junctions in read alignments. Ideally, we recommend to use the same genome seqeunce file that was used for alignment to be used for bambu. ```{r} # get path to fasta file genomeSequenceData <- query(ExperimentHub(), c("NanoporeRNA", "GRCh38","FASTA")) genomeSequence <- genomeSequenceData[["EH7260"]] ``` As an option, users can also choose to use a BSgenome object: ```{r} library(BSgenome.Hsapiens.NCBI.GRCh38) genomeSequenceBsgenome <- BSgenome.Hsapiens.NCBI.GRCh38 ``` ### Genome annotations (bambu annotations object/ gtf file / TxDb object) *bambu* also requires a reference transcript annotations object, which is used to correct read alignments, to train a model for transcript discovery, and for quantification. The annotation object can be created from a gtf file: ```{r} gtf.file <- system.file("extdata", "Homo_sapiens.GRCh38.91_chr9_1_1000000.gtf", package = "bambu") annotation <- prepareAnnotations(gtf.file) ``` The annotation object can also be created from a TxDb object: ```{r} txdb <- system.file("extdata", "Homo_sapiens.GRCh38.91_chr9_1_1000000.gtf", package = "bambu") annotation <- prepareAnnotations(txdb) ``` The annotation object can be stored and used again for re-running bambu. Here we will used the annotation object from the *NanoporeRNASeq* package that wasis prepared from a gtf file using the \code{prepareAnnotations} function in by function in *bambu*. ```{r} data("HsChr22BambuAnnotation") HsChr22BambuAnnotation ``` ## Transcript discovery and quantification ### Running bambu Next we apply *bambu* on the input data (bam files, annotations, genomeSequence). Bambu will perform isoform discovery to extend the provided annotations, and then quantify the transcript expression from these extended annotations using an Expectation-Maximisation algorithm. Here we will use 1 core, which can be changed to process multiple files in parallel. ```{r, results = "hide"} se <- bambu(reads = bamFile, annotations = HsChr22BambuAnnotation, genome = genomeSequence, ncore = 1) ``` ```{r} se ``` Optionally, users can choose to apply *bambu* to do transcript discovery only ```{r, results = "hide"} seNoquant <- bambu(reads = bamFile, annotations = HsChr22BambuAnnotation, genome = genomeSequence, quant = FALSE) ``` Optionally, users can choose to apply *bambu* to do quantification only (without isoform discovery) ```{r, results = "hide"} seUnextended <- bambu(reads = bamFile, annotations = HsChr22BambuAnnotation, genome = genomeSequence, discovery = FALSE) ``` ```{r} seUnextended ``` ### Running multiple samples {#running-multiple-samples} If you have multiple replicates for a sample, or plan to do comparative analysis between conditions, it may be beneficial to run all samples together instead of individually. This can be done by providing a vector of paths to all the bam files you want to analyze together. ```{r, results = "hide"} bamFiles <- Rsamtools::BamFileList(NanoporeData[["EH3808"]], NanoporeData[["EH3809"]],NanoporeData[["EH3810"]], NanoporeData[["EH3811"]], NanoporeData[["EH3812"]], NanoporeData[["EH3813"]]) se.multiSample <- bambu(reads = bamFiles, annotations = HsChr22BambuAnnotation, genome = genomeSequence) ``` The advantage of running samples together include: novel transcripts that are identified in multiple samples are assigned unified IDs, enabling comparative analysis between different samples. This is especially important for downstream differential expression analysis when looking at novel transcripts. Running multiple samples can be multithreaded. While running multiple samples, by default, *bambu* will train a model separately on each sample and score novel transcripts in each sample separately. If you need to combine samples in multiple configurations (for example different pairwise comparisons) we would recommend using the intermediate rcFiles to save processing time. ### Modulating the sensitivity of discovery (pre and post analysis) {#modulating-the-sensitivity} When doing transcript discovery there is a balance between sensitivity (the number of real novel transcripts that are detected) and the precision (how many of the novel transcripts are real). To control this balance, *bambu* uses the novel discovery rate (NDR) as the main parameter. The NDR threshold approximates the proportion of novel candidates output by bambu, relative to the number of known transcripts it found, i.e., an NDR of 0.1 would mean that 10% of all transcripts passing the threshold are classified as novel. If you are using a genome where you expect a high amount of novel transcripts, a higher NDR is recommended so that these novel transcripts are not missed. Conversely if you are using a well annotated genome, we recommend a lower NDR threshold to reduce the presence of false positives. By default the NDR threshold is automatically chosen for the user based on predicted level of annotation completeness when compared to the default model trained on human reference annotations (Hg38). To manually select an NDR value, use the NDR argument in *bambu*: ```{r, results = "hide"} se.NDR_0.3 <- bambu(reads = bamFiles, annotations = HsChr22BambuAnnotation, genome = genomeSequence, NDR = 0.3) ``` Alternatively transcript discovery can be run without thresholds, producing a GRangesList annotation object with all transcripts scored with its NDR score. Note that this means turning quant = FALSE in running *bambu*. The annotations can be filtered by their NDR score (see example below), read count and gene read proportion between the discovery and quantification steps or used for other types of analysis. ```{r, results = 'hide'} newAnnotations <- bambu(reads = bamFiles, annotations = HsChr22BambuAnnotation, genome = genomeSequence, NDR = 1, quant = FALSE) annotations.filtered <- newAnnotations[(!is.na(mcols(newAnnotations)$NDR) & mcols(newAnnotations)$NDR<0.1) | is.na(mcols(newAnnotations)$NDR)] se.NDR_1 <- bambu(reads = bamFiles, annotations = annotations.filtered, genome = genomeSequence, NDR = 1, discovery = FALSE) ``` Additionally there are other thresholds that advanced users can access through opt.discovery when running *bambu* (see arguments). ### Visualise results {#visualise-results} *bambu* provides functions to visualise and explore the results. When multiple samples are used, we can visualise the correlation and clustering of all samples with a heatmap: ```{r, fig.width = 8, fig.height = 6} library(ggplot2) plotBambu(se.multiSample, type = "heatmap") ``` Additionally, we can also visualise the correlation with a 2-dimmensional PCA plot. ```{r, fig.width = 8, fig.height = 6} plotBambu(se.multiSample, type = "pca") ``` In addition to visualising the correlation between samples, *bambu* also provide a function to visualise the extended annotation and expression estimation for individual genes. Here we look at gene ENSG00000099968 and visualise the transcript coordinates for annotated and novel isoforms and expression levels for these isoforms across all samples. ```{r, fig.width = 8, fig.height = 10} plotBambu(se.multiSample, type = "annotation", gene_id = "ENSG00000099968") ``` ### Obtain gene expression estimates from transcript expression {#gene-expression} For the downstream analysis, we will add the condition of interest to the \code{colData} object that describes the samples. Here we are interested in a comparison of the 2 cell lines: ```{r} colData(se.multiSample)$condition <- as.factor(SGNexSamples$cellLine) ``` To obtain the accurate gene expression estimates which uses all reads that can be assigned to each gene (including reads that are incompatible with all existing annotations) you can run the following command. Looking at the output, we can see there are novel genes identified as well ```{r} seGene.multiSample <- transcriptToGeneExpression(se.multiSample) seGene.multiSample ``` We can again use the \code{plotBambu} function to visualise the gene expression data across the 6 samples with a heatmap or PCA plot. As expected, samples from the same cell line showed higher correlation than across the cell lines. ```{r, fig.width = 8, fig.height = 6} colData(seGene.multiSample)$groupVar <- SGNexSamples$cellLine plotBambu(seGene.multiSample, type = "heatmap") ``` ### Save data (gtf/text) {#save-data} *bambu* includes a function to write the extended annotations, the transcript and the gene expression estimates that include any newly discovered genes and transcripts to text files. ```{r} save.dir <- tempdir() writeBambuOutput(se.multiSample, path = save.dir, prefix = "NanoporeRNASeq_") ``` *bambu* also includes a function that only exports the extended annotations to gtf file: ```{r} save.file <- tempfile(fileext = ".gtf") writeToGTF(rowRanges(se.multiSample), file = save.file) ``` ## Advanced usages for different use cases {#advanced-usages} ### Using a pre-trained model {#use-pre-trained-model} By default, *bambu* requires at least 1000 transcripts from the annotations to be detected in a sample in order to train a sample specific model. In use cases where this is not possible *bambu* will instead use a default pre-trained model to calculate the transcript probability score (TPS) for each read class. Users can force this behavior if they believe their annotations are not sufficient for sample specific training (for example if they suspect a high proportion of real novel transcripts are present in their sample). This is advantageous when you want NDR calibration without the impacts of a model trained using low quality annotations. ```{r, results = 'hide'} se <- bambu(reads = bamFiles, annotations = HsChr22BambuAnnotation, genome = genomeSequence, opt.discovery = list(fitReadClassModel = FALSE)) ``` The default pre-trained model was trained on SGNex_HepG2_directRNA_replicate5_run1 and has the following characteristics: Genome: Homo_sapiens.GRCh38.dna_sm.primary_assembly
Annotations: Homo_sapiens.GRCh38.91
Read count: 7,861,846
Technology: Nanopore (ONT)
Library preparation: directRNA
Base Calling Accuracy: 79%
Average Read Length: 1093
We have found the pre-trained model works successfully across species borders (on Arabidopsis thaliana) and on different technologies (PacBio), with only small decreases in performance compared to using a sample specific model. The pre-trained model is not always effective in samples with large differences in sequencing quality or if the library preparation results in biases in the overall structure of the transcriptome. In this case, we would recommend training a new model using similar data from a different sample that has quality reference annotations (See [Training a model on another species/dataset and applying it](#train-model)). For the example in this vignette, the data used is rather small, hence the pre-trained model is used by default. So you won't see differences in the output, but you may try it out with more data (see [SG-NEx data](https://github.com/GoekeLab/sg-nex-data)) ### De-novo transcript discovery {#de-novo} In cases where the organism does not yet have reference annotations, or unreliable annotations, *bambu* can be run in de-novo mode. In de-novo mode, *bambu* does not train a model, and instead uses the pre-trained model to classify novel transcripts (see [Using a pre-trained model](#use-pretrained-model). To learn how to train a new model for a more closely related organism/sample see [Training a model on another species/dataset and applying it](#train-model). Without annotations *bambu* is unable to calibrate the NDR output, nor be able to recommend a threshold and will instead use the TPS as the thresholded value. Therefore you should supply a manual NDR threshold ([Modulating the sensitivity of discovery (pre and post analysis)](#modulating-the-sensitivity)) and note that the precision of the output is unlikely to linearly match an applied threshold. The TPS threshold used is (> 1-NDR). If an NDR is not provided, a default NDR threshold of <0.1 is used (an effective TPS threshold of > 0.9). As in [Modulating the sensitivity of discovery (pre and post analysis)](#modulating-the-sensitivity) an NDR of 1 can be provided to output all possible read classes with their TPS scores ```{r, results = 'hide'} novelAnnotations <- bambu(reads = bamFiles, annotations = NULL, genome = genomeSequence, NDR = 1, quant = FALSE) ``` ### Storing and using preprocessed files (rcFiles) The first step of *bambu* involves the construction of read classes which is a large fraction of the running time. This could be time-consuming if we want to perform transcript discovery & quantification multiple times on the same dataset using different configurations (eg. NDR, or combinations of samples), especially when the sample is large. To mitigate this, we can store the read class information as read class files (rcFiles) during a *bambu* run. Then they can be used as an input argument in the *bambu* main function for the subsequent *bambu* run. ```{r, eval = FALSE} #se <- bambu(rcFile = rcFiles, annotations = annotations, genome = fa.file) ``` rcFiles can be generated in two ways, either as a direct output of the bambu() function when quant and discovery are FALSE, or as written outputs when a path is provided to the rcOutdir argument. When rcFiles are output using rcOutdir this is done using BiocFileCache. For more details on how to access, use, and identify these files see [here](https://bioconductor.org/packages/release/bioc/html/BiocFileCache.html). A short example is shown below. Example using rcOutDir to produce preprocessed files ```{r, results = 'hide'} se <- bambu(reads = bamFiles, rcOutDir = save.dir, annotations = HsChr22BambuAnnotation, genome = genomeSequence) ``` This will store a preprocessed rcFile in the provided directory for each sample file provided to reads. To access these files for future use, we recommend using the BioCFileCache package which provides the metadata needed to identify the sample. ```{r} library(BiocFileCache) bfc <- BiocFileCache(save.dir, ask = FALSE) info <- bfcinfo(bfc) ``` The info object is a tibble which associates the filename (fpath) with the sample (rname) to help you identify which .rds file you need. ```{r} info # running bambu using the first file se <- bambu(reads = info$rpath[1:3], annotations = HsChr22BambuAnnotation, genome = genomeSequence) ``` This output is also generated when both quant and discovery are set to false in a list form indexed by sample. ```{r} se <- bambu(reads = bamFiles, annotations = HsChr22BambuAnnotation, genome = genomeSequence, discovery = FALSE, quant = FALSE) ``` As this is an intermediate object it will not be suitable to use for general use cases. We will document the object below for any potential advanced use cases that may arise. ```{r} rowData(se[[1]]) ``` ### Tracking read-to-transcript assignment {#track-read-assignment} Some use cases require knowing which individual reads support specific transcripts (novel and annotated). By default this feature is off due to the memory overhead it introduces but can be turned on using the trackReads argument. The output has three columns: read_id, a list of indices of equal matches, a list of indices of compatible matches. These indices match the annotations found in rowRanges(se) ```{r} se <- bambu(reads = bamFiles, annotations = HsChr22BambuAnnotation, genome = genomeSequence, trackReads = TRUE) metadata(se)$readToTranscriptMaps[[1]] ``` ### Training a model on another species/dataset and applying it {#train-model} In situations where training is not or cannot be performed, and the default model is also not suitable for the sample (the sample is generated from a different technology, species, condition, etc), *bambu* provides the option to train a new model, if well annotated similar data is available. For example one might train a model on arabidopsis to apply to an unannotated plant sample. ```{r, eval = FALSE} # first train the model using a related annotated dataset from .bam se = bambu(reads = sample1.bam, annotations = annotations, genome = fa.file, discovery = FALSE, quant = FALSE, opt.discovery = list(returnModel = TRUE)) # note that discovery and quant need to be set to FALSE, alternatively you can have them set to TRUE and retrieve the model from the rcFile as long as returnModel = TRUE. # newDefaultModel = metadata(se[[1]])$model # [[1]] will select the model trained on the first sample #alternatively train the model using an rcFile rcFile <- readRDS(pathToRcFile) newDefaultModel = trainBambu(rcFile) # use the trained model on another sample # sample2.bam and fa.file2 represent the aligned reads and genome for the poorly annotated sample se <- bambu(reads = sample2.bam, annotations = NULL, genome = fa.file2, opt.discovery = list(defaultModels = newDefaultModel, fitReadClassModel = FALSE)) #trainBambu Arguments rcFile <- NULL, min.readCount = 2, nrounds = 50, NDR.threshold = 0.1, verbose = TRUE ``` ### Including single exons {#include-single-exons} By default *bambu* does not report single exon transcripts because they are known to have a high frequency of false positives and do not have splice junctions that are used by *bambu* to distinguish read classes. Nevertheless *bambu* trains a separate model on single-exon transcripts, and these predictions can be accessed and included in the annotations. ```{r} se <- bambu(reads = bamFiles, annotations = HsChr22BambuAnnotation, genome = genomeSequence, opt.discovery = list(min.txScore.singleExon = 0)) ``` ## Downstream analysis ### Identifying differentially expressed genes {#DESeq2} One of the most common tasks when analysing RNA-Seq data is the analysis of differential gene expression across a condition of intertest. Here we use *DESeq2* to find the differentially expressed genes between MCF7 and K562 cell lines. Similar to using results from *Salmon*, estimates from *bambu* will first be rounded. ```{r} library(DESeq2) dds <- DESeqDataSetFromMatrix(round(assays(seGene.multiSample)$counts), colData = colData(se.multiSample), design = ~ condition) dds.deseq <- DESeq(dds) deGeneRes <- DESeq2::results(dds.deseq, independentFiltering = FALSE) head(deGeneRes[order(deGeneRes$padj),]) ``` A quick summary of differentially expressed genes ```{r} summary(deGeneRes) ``` We can also visualise the MA-plot for differentially used isoforms using \code{plotMA(deGeneRes)}. However, visualizing the MA-plots using the original log-fold change results will be affected by the noise associated with log2 fold changes from low count genes without requiring arbitrary filtering thresholds. As recommended in the [*DESeq2* tutorial](http://bioconductor.org/packages/devel/bioc/vignettes/DESeq2/inst/doc/DESeq2.html#alternative-shrinkage-estimators). we applied the same shrinkage to effect sizes to improve the visualization. ```{r, fig.width = 8, fig.height = 6} library(apeglm) resLFC <- lfcShrink(dds.deseq, coef = "condition_MCF7_vs_K562", type = "apeglm") plotMA(resLFC, ylim = c(-3,3)) ``` ### Identifying differential transcript usage {#DEXSeq} We used *DEXSeq* to detect alternative used isoforms. ```{r} library(DEXSeq) dxd <- DEXSeqDataSet(countData = round(assays(se.multiSample)$counts), sampleData = as.data.frame(colData(se.multiSample)), design = ~sample + exon + condition:exon, featureID = rowData(se.multiSample)$TXNAME, groupID = rowData(se.multiSample)$GENEID) dxr <- DEXSeq(dxd) head(dxr) ``` We can visualize the MA-plot ```{r,fig.width = 8, fig.height = 6} plotMA(dxr, cex = 0.8 ) ``` # Getting help {#get-help} Questions and issues can be raised at the Bioconductor support site (once bambu is available through bioconductor): https://support.bioconductor.org. Please tag your your posts with bambu. Alternatively, issues can be raised at the bambu Github repository:https://github.com/GoekeLab/bambu. # Citing bambu {#cite-bambu} A manuscript describing bambu is currently in preparation. If you use bambu for your research, please cite using the following doi: 10.5281/zenodo.3900025. # Session Information {#session-info} ```{r} sessionInfo() ```