gettRNAstructureSeqs: tRNA structures and sequences

gettRNAstructureGRangesR Documentation

tRNA structures and sequences

Description

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"]].

Usage

gettRNAstructureGRanges(x, structure = "")

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

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

## S4 method for signature '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")
gettRNAstructureSeqs(gr, structure = "anticodonLoop")
gettRNABasePairing(gr[1:10])

FelixErnst/tRNA documentation built on April 1, 2024, 2:39 p.m.