predictNMD | R Documentation |
Predict NMD sensitivity on mRNA transcripts
predictNMD(x, ..., cds = NULL, NMD_threshold = 50, progress_bar = TRUE)
x |
Can be a GRanges object containing exon and CDS transcript features in GTF format. Can be a GRangesList object containing exon features for a list of transcripts.If so, 'cds' argument have to be provided. Can be a GRanges object containing exon features for a transcript. If so, 'cds' argument have to be provided. |
... |
Logical conditions to pass to dplyr::filter to subset transcripts for analysis. Variables are metadata information found in 'x' and multiple conditions can be provided delimited by comma. Example: transcript_id == "transcript1" |
cds |
If 'x' is a GRangesList object, 'cds' has to be a GRangesList containing CDS features for the list of transcripts in 'x'. List names in 'x' and 'cds' have to match. If 'x' is a GRanges object, 'cds' has to be a GRanges containing CDS features for the transcript in 'x'. |
NMD_threshold |
Minimum distance of stop_codon to last exon junction (EJ) which triggers NMD. Default = 50bp |
progress_bar |
Whether to display progress Default = TRUE |
Dataframe with prediction of NMD sensitivity and NMD features:
is_NMD: logical value in prediciting transcript sensitivity to NMD
stop_to_lastEJ: Integer value of the number of bases between the first base of the stop_codon to the last base of EJ. A positive value indicates that the last EJ is downstream of the stop_codon.
num_of_down_EJs: Number of EJs downstream of the stop_codon.
‘3_UTR_length': Length of 3’ UTR
Fursham Hamid
## ---------------------------------------------------------------------
## EXAMPLE USING SAMPLE DATASET
## ---------------------------------------------------------------------
# Load datasets
data(new_query_gtf, query_exons, query_cds)
## Using GTF GRanges as input
predictNMD(new_query_gtf)
### Transcripts for analysis can be subsetted using logical conditions
predictNMD(new_query_gtf, transcript_id == "transcript1")
predictNMD(new_query_gtf,
transcript_id %in% c("transcript1", "transcript3"))
## Using exon and CDS GRangesLists as input
predictNMD(query_exons, cds = query_cds)
predictNMD(query_exons, cds = query_cds, transcript_id == "transcript3")
## Using exon and CDS GRanges as input
predictNMD(query_exons[[3]], cds = query_cds[[3]])
## ---------------------------------------------------------------------
## EXAMPLE USING TRANSCRIPT ANNOTATION
## ---------------------------------------------------------------------
library(AnnotationHub)
## Retrieve GRCm38 trancript annotation
ah <- AnnotationHub()
GRCm38_gtf <- ah[["AH60127"]]
## Run tool on specific gene family
predictNMD(GRCm38_gtf, gene_name == "Ptbp1")
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.