genomeToProtein | R Documentation |
Map positions along the genome to positions within the protein sequence if
a protein is encoded at the location. The provided coordinates have to be
completely within the genomic position of an exon of a protein coding
transcript (see genomeToTranscript()
for details). Also, the provided
positions have to be within the genomic region encoding the CDS of a
transcript (excluding its stop codon; soo transcriptToProtein()
for
details).
For genomic positions for which the mapping failed an IRanges
with
negative coordinates (i.e. a start position of -1) is returned.
genomeToProtein(x, db, proteins = NA, exons = NA, transcripts = NA)
x |
|
db |
|
proteins |
|
exons |
|
transcripts |
|
genomeToProtein
combines calls to genomeToTranscript()
and
transcriptToProtein()
.
An IRangesList
with each element representing the mapping of one of the
GRanges
in x
(i.e. the length of the IRangesList
is length(x)
).
Each element in IRanges
provides the coordinates within the protein
sequence, names being the (Ensembl) IDs of the protein. The ID of the
transcript encoding the protein, the ID of the exon within which the
genomic coordinates are located and its rank in the transcript are provided
in metadata columns "tx_id"
, "exon_id"
and "exon_rank"
. Metadata
columns "cds_ok"
indicates whether the length of the CDS matches the
length of the encoded protein. Coordinates for which cds_ok = FALSE
should
be taken with caution, as they might not be correct. Metadata columns
"seq_start"
, "seq_end"
, "seq_name"
and "seq_strand"
provide the
provided genomic coordinates.
For genomic coordinates that can not be mapped to within-protein sequences
an IRanges
with a start coordinate of -1 is returned.
Johannes Rainer
Other coordinate mapping functions:
cdsToTranscript()
,
genomeToTranscript()
,
proteinToGenome()
,
proteinToTranscript()
,
transcriptToCds()
,
transcriptToGenome()
,
transcriptToProtein()
library(EnsDb.Hsapiens.v86)
## Restrict all further queries to chromosome x to speed up the examples
edbx <- filter(EnsDb.Hsapiens.v86, filter = ~ seq_name == "X")
## In the example below we define 4 genomic regions:
## 630898: corresponds to the first nt of the CDS of ENST00000381578
## 644636: last nt of the CDS of ENST00000381578
## 644633: last nt before the stop codon in ENST00000381578
## 634829: position within an intron.
gnm <- GRanges("X", IRanges(start = c(630898, 644636, 644633, 634829),
width = c(5, 1, 1, 3)))
res <- genomeToProtein(gnm, edbx)
## The result is an IRangesList with the same length as gnm
length(res)
length(gnm)
## The first element represents the mapping for the first GRanges:
## the coordinate is mapped to the first amino acid of the protein(s).
## The genomic coordinates can be mapped to several transcripts (and hence
## proteins).
res[[1]]
## The stop codon is not translated, thus the mapping for the second
## GRanges fails
res[[2]]
## The 3rd GRanges is mapped to the last amino acid.
res[[3]]
## Mapping of intronic positions fail
res[[4]]
## Meanwhile, this function can be called in parallel processes if you preload
## the protein, exons and transcripts database.
proteins <- proteins(edbx)
exons <- exonsBy(edbx)
transcripts <- transcripts(edbx)
genomeToProtein(gnm, edbx, proteins = proteins, exons = exons, transcripts = transcripts)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.