GRanges
Example of importing a tRNAscan-SE output for sacCer3 as a GRanges object
tRNAscanImport 1.27.0
tRNAscan-SE (Lowe and Eddy 1997) can be used for prediction of tRNA genes in whole
genomes based on sequence context and calculated structural features. Many tRNA
annotations in genomes contain or are based on information generated by
tRNAscan-SE, for example the current SGD reference genome sacCer3 for
Saccharomyces cerevisiae. However, not all available information from
tRNAscan-SE end up in the genome annotation. Among these are for example
structural information, additional scores and the information, whether the
conserved CCA-end is encoded in the genomic DNA. To work with this complete set
of information, the tRNAscan-SE output can be parsed into a more accessible
GRanges object using tRNAscanImport
.
The default tRNAscan-SE output, either from running tRNAscan-SE (Lowe and Eddy 1997) locally or retrieving the output from the gtRNADb (Chan and Lowe 2016), consist of a formatted text document containing individual text blocks per tRNA delimited by an empty line.
library(tRNAscanImport)
yeast_file <- system.file("extdata",
file = "yeast.tRNAscan",
package = "tRNAscanImport")
# output for sacCer3
# Before
readLines(con = yeast_file, n = 7L)
## [1] "chrI.trna1 (139152-139254)\tLength: 103 bp"
## [2] "Type: Pro\tAnticodon: TGG at 33-35 (139184-139186)\tScore: 62.1"
## [3] "Possible intron: 37-67 (139188-139218)"
## [4] "HMM Sc=37.90\tSec struct Sc=24.20"
## [5] " * | * | * | * | * | * | * | * | * | * | "
## [6] "Seq: GGGCGTGTGGTCTAGTGGTATGATTCTCGCTTTGGGcgacttcctgattaaacaggaagacaaagcaTGCGAGAGGcCCTGGGTTCAATTCCCAGCTCGCCCC"
## [7] "Str: >>>>>.>..>>>.........<<<.>>>>>......................................<<<<<.....>>>>>.......<<<<<<.<<<<<."
GRanges
To access the information in a BioC context the import as a GRanges object
comes to mind. import.tRNAscanAsGRanges()
performs this task by evaluating
each text block using regular expressions.
# output for sacCer3
# After
gr <- import.tRNAscanAsGRanges(yeast_file)
head(gr, 2)
## GRanges object with 2 ranges and 18 metadata columns:
## seqnames ranges strand | no tRNA_length tRNA_type
## <Rle> <IRanges> <Rle> | <integer> <integer> <character>
## [1] chrI 139152-139254 + | 1 72 Pro
## [2] chrI 166267-166339 + | 2 73 Ala
## tRNA_anticodon tRNA_anticodon.start tRNA_anticodon.end tRNAscan_score
## <character> <integer> <integer> <numeric>
## [1] TGG 33 35 62.1
## [2] TGC 34 36 76.0
## tRNA_seq tRNA_str tRNA_CCA.end
## <DNAStringSet> <DotBracketStringSet> <logical>
## [1] GGGCGTGTGG...AGCTCGCCCC <<<<<.<..<...>>>.>>>>>. FALSE
## [2] GGGCACATGG...GTTGCGTCCA <<<<.<<..<...>>>>.>>>>. FALSE
## tRNAscan_potential.pseudogene tRNAscan_intron.start tRNAscan_intron.end
## <logical> <integer> <integer>
## [1] FALSE 139188 139218
## [2] FALSE <NA> <NA>
## tRNAscan_intron.locstart tRNAscan_intron.locend tRNAscan_hmm.score
## <integer> <integer> <numeric>
## [1] 37 67 37.9
## [2] <NA> <NA> 53.4
## tRNAscan_sec.str.score tRNAscan_infernal
## <numeric> <numeric>
## [1] 24.2 NA
## [2] 22.6 NA
## -------
## seqinfo: 17 sequences from an unspecified genome; no seqlengths
# Any GRanges passing this, can be used for subsequent function
istRNAscanGRanges(gr)
## [1] TRUE
The result can be used directly in R or saved as gff3/fasta file for further use, including processing the sequences for HTS read mapping or statistical analysis on tRNA content of the analyzed genome.
library(Biostrings)
library(rtracklayer)
# suppressMessages(library(rtracklayer, quietly = TRUE))
# Save tRNA sequences
writeXStringSet(gr$tRNA_seq, filepath = tempfile())
# to be GFF3 compliant use tRNAscan2GFF
gff <- tRNAscan2GFF(gr)
export.gff3(gff, con = tempfile())
The tRNAscan-SE information can be visualized using the gettRNAFeaturePlots()
function of the tRNA
package, returning a named list of ggplot2 plots, which
can be plotted or further modified.
Alternatively, gettRNASummary()
returns the aggregated information for
further use.
# tRNAscan-SE output for hg38
human_file <- system.file("extdata",
file = "human.tRNAscan",
package = "tRNAscanImport")
# tRNAscan-SE output for E. coli MG1655
eco_file <- system.file("extdata",
file = "ecoli.tRNAscan",
package = "tRNAscanImport")
# import tRNAscan-SE files
gr_human <- import.tRNAscanAsGRanges(human_file)
gr_eco <- import.tRNAscanAsGRanges(eco_file)
# get summary plots
grl <- GRangesList(Sce = gr,
Hsa = gr_human,
Eco = gr_eco)
plots <- gettRNAFeaturePlots(grl)
plots$length