gettRNAstructureGRanges returns a list of GRanges describing the boundaries of tRNA structures as extracted from the dot bracket annotation. The dot bracket annotation is parsed using gettRNABasePairing, which internally uses getBasePairing.

gettRNAstructureSeq returns split or partial tRNA sequences based on the structure information. Variations in the ength of structure features can be padded to retrieve sequences of equal length. If sequences are joined by setting joinCompletely = FALSE, the boundaries of the tRNA structure are stored in the result as metadata. They can be accessesed as an IRanges object by using metadata(seq)[["tRNA_structures"]].

gettRNAstructureGRanges(x, structure = "")

gettRNAstructureSeqs(
  x,
  structure = "",
  joinCompletely = FALSE,
  joinFeatures = FALSE,
  padSequences = TRUE
)

# S4 method for GRanges
gettRNAstructureSeqs(
  x,
  structure = "",
  joinCompletely = FALSE,
  joinFeatures = FALSE,
  padSequences = TRUE
)

# S4 method for GRanges
gettRNAstructureGRanges(x, structure = "")

Arguments

x

a GRanges object with tRNA information. It has to pass the istRNAGRanges function.

structure

optional parameter for returning just partial sequences. The following values are accepted: anticodonStem, Dprime5, DStem, Dloop, Dprime3, acceptorStem, anticodonloop, variableLoop, TStem, Tloop, discriminator. (default: structure = "")

joinCompletely

Should the sequence parts, which are to be returned, be joined into one sequence? (default: joinCompletely = FALSE)) Setting this to TRUE excludes joinFeatures be set to TRUE as well. In addition, joinCompletely = TRUE uses automatically all sequence structures.

joinFeatures

Should the sequence parts, which are to be returned and are from the same structure type, be joined into one sequence? (default: joinCompletely = FALSE)) Setting this to TRUE excludes joinCompletely be set to TRUE as well. joinCompletely takes precedence.

padSequences

parameter whether sequences of the same type should be returned with the same length. For stems missing positions will be filled up in the middle, for loops at the ends. (default: padSequences = TRUE). If joinCompletely == TRUE this is set to TRUE automatically.

Value

a list of GRanges or DNAStringSet objects. In case joinCompletly is set to TRUE a single DNAStringSet is returned.

Examples

data("gr", package = "tRNA")
gettRNAstructureGRanges(gr, structure = "anticodonLoop")
#> $anticodonLoop
#> IRanges object with 299 ranges and 0 metadata columns:
#>           start       end     width
#>       <integer> <integer> <integer>
#>   TGG        31        37         7
#>   TGC        32        38         7
#>   CAA        31        37         7
#>   AGA        31        37         7
#>   TAA        31        37         7
#>   ...       ...       ...       ...
#>   CAT        32        38         7
#>   GAA        31        37         7
#>   TTA        31        37         7
#>   TAC        32        38         7
#>   CAT        32        38         7
#> 
gettRNAstructureSeqs(gr, structure = "anticodonLoop")
#> $anticodonLoop
#> RNAStringSet object of length 299:
#>       width seq                                             names               
#>   [1]     7 UUUGGGU                                         TGG
#>   [2]     7 CUUGCAA                                         TGC
#>   [3]     7 UUCAAGC                                         CAA
#>   [4]     7 UUAGAAA                                         AGA
#>   [5]     7 CUUAAGA                                         TAA
#>   ...   ... ...
#> [295]     7 CUCAUAA                                         CAT
#> [296]     7 UUGAAGA                                         GAA
#> [297]     7 UUUUAGU                                         TTA
#> [298]     7 UUUACAC                                         TAC
#> [299]     7 GUCAUGA                                         CAT
#> 
gettRNABasePairing(gr[1:10])
#> SplitDotBracketDataFrameList of length 10
#> [[1]]
#> DotBracketDataFrame with 71 rows and 4 columns
#>           pos   forward   reverse   character
#>     <integer> <integer> <integer> <character>
#> 1           1         1        70           <
#> 2           2         2        69           <
#> 3           3         3        68           <
#> 4           4         4        67           <
#> 5           5         5        66           <
#> ...       ...       ...       ...         ...
#> 67         67        67         4           >
#> 68         68        68         3           >
#> 69         69        69         2           >
#> 70         70        70         1           >
#> 71         71         0         0           .
#> 
#> ...
#> <9 more elements>