findMapORFs {ORFik} | R Documentation |
Finds ORFs on the sequences of interest, but returns relative positions to the positions of 'grl' argument. For example, 'grl' can be exons of known transcripts (with genomic coordinates), and 'seq' sequences of those transcripts, in that case, [findMapORFs()] will return genomic coordinates of ORFs found on transcript sequences.
findMapORFs(grl, seqs, startCodon = startDefinition(1), stopCodon = stopDefinition(1), longestORF = FALSE, minimumLength = 0)
grl |
( |
seqs |
(DNAStringSet or character) DNA sequences to search for Open Reading Frames. |
startCodon |
(character) Possible START codons to search for. Check [startDefinition()] for helper function. |
stopCodon |
(character) Possible STOP codons to search for. Check [stopDefinition()] for helper function. |
longestORF |
(logical) Default FALSE. When TRUE will only report ORFs that are longest, all smaller overlapping ORFs will be ignored. When FALSE will report all possible ORFs in all three reading frames. |
minimumLength |
(integer) Default is 0. Minimum length of ORF, without counting 3bp for START and STOP codons. For example minimumLength = 8 will result in size of ORFs to be at least START + 8*3 (bp) + STOP. Use this param to restrict search. |
This function assumes that 'seq' is in widths relative to 'grl', and that their orders match.
A GRangesList of ORFs.
[findORFs()], [findORFsFasta()], [startDefinition()], [stopDefinition()]
Other findORFs: findORFsFasta
,
findORFs
, startDefinition
,
stopDefinition
# This sequence has ORFs at 1-9 and 4-9 seqs <- c("ATGATGTAA") # the dna sequence findORFs(seqs) # lets assume that this sequence comes from two exons as follows gr <- GRanges(seqnames = rep("1", 2), # chromosome 1 ranges = IRanges(start = c(21, 10), end = c(23, 15)), strand = rep("-", 2), names = rep("tx1", 2)) grl <- GRangesList(tx1 = gr) findMapORFs(grl, seqs) # ORFs are properly mapped to its genomic coordinates grl <- c(grl, grl) names(grl) <- c("tx1", "tx2") findMapORFs(grl, c(seqs, seqs))