Title: | Mapping Alternative Splicing Events to pRoteins |
---|---|
Description: | This package provides functionalities for downstream analysis, annotation and visualizaton of alternative splicing events generated by rMATS. |
Authors: | Diogo F.T. Veiga [aut, cre] |
Maintainer: | Diogo F.T. Veiga <[email protected]> |
License: | MIT + file LICENSE |
Version: | 1.25.0 |
Built: | 2024-11-19 03:59:36 UTC |
Source: | https://github.com/bioc/maser |
Retrieve annotation of splicing events from a maser object.
## S4 method for signature 'Maser' annotation(object, type = c("A3SS", "A5SS", "SE", "RI", "MXE"))
## S4 method for signature 'Maser' annotation(object, type = c("A3SS", "A5SS", "SE", "RI", "MXE"))
object |
a maser object. |
type |
a character indicating the splice type. Possible values
are |
a data.frame.
path <- system.file("extdata", file.path("MATS_output"), package = "maser") hypoxia <- maser(path, c("Hypoxia 0h", "Hypoxia 24h")) head(annotation(hypoxia, "SE"))
path <- system.file("extdata", file.path("MATS_output"), package = "maser") hypoxia <- maser(path, c("Hypoxia 0h", "Hypoxia 24h")) head(annotation(hypoxia, "SE"))
Query available human protein features in UniprotKB.
availableFeaturesUniprotKB()
availableFeaturesUniprotKB()
a data.frame.
head(availableFeaturesUniprotKB(), 10)
head(availableFeaturesUniprotKB(), 10)
Boxplots of PSI distributions by splicing type.
boxplot_PSI_levels(events, type = c("A3SS", "A5SS", "SE", "RI", "MXE"))
boxplot_PSI_levels(events, type = c("A3SS", "A5SS", "SE", "RI", "MXE"))
events |
a maser object. |
type |
character indicating splice type. Possible values are
|
a ggplot object.
path <- system.file("extdata", file.path("MATS_output"), package = "maser") hypoxia <- maser(path, c("Hypoxia 0h", "Hypoxia 24h")) hypoxia_filt <- filterByCoverage(hypoxia, avg_reads = 5) boxplot_PSI_levels(hypoxia_filt, type = "RI")
path <- system.file("extdata", file.path("MATS_output"), package = "maser") hypoxia <- maser(path, c("Hypoxia 0h", "Hypoxia 24h")) hypoxia_filt <- filterByCoverage(hypoxia, avg_reads = 5) boxplot_PSI_levels(hypoxia_filt, type = "RI")
Retrieve raw read counts values from a maser object.
## S4 method for signature 'Maser' counts(object, type = c("A3SS", "A5SS", "SE", "RI", "MXE"))
## S4 method for signature 'Maser' counts(object, type = c("A3SS", "A5SS", "SE", "RI", "MXE"))
object |
a maser object. |
type |
a character indicating the splice type. Possible values
are |
a matrix.
path <- system.file("extdata", file.path("MATS_output"), package = "maser") hypoxia <- maser(path, c("Hypoxia 0h", "Hypoxia 24h")) head(counts(hypoxia, "SE"))
path <- system.file("extdata", file.path("MATS_output"), package = "maser") hypoxia <- maser(path, c("Hypoxia 0h", "Hypoxia 24h")) head(counts(hypoxia, "SE"))
Visualization of splicing events annotation using an interactive data table.
display(events, type = c("A3SS", "A5SS", "SE", "RI", "MXE"))
display(events, type = c("A3SS", "A5SS", "SE", "RI", "MXE"))
events |
a maser object. |
type |
character indicating splice type. Possible values are
|
a datatables object.
path <- system.file("extdata", file.path("MATS_output"), package = "maser") hypoxia <- maser(path, c("Hypoxia 0h", "Hypoxia 24h")) hypoxia_filt <- filterByCoverage(hypoxia, avg_reads = 5) hypoxia_top <- topEvents(hypoxia_filt) display(hypoxia_top, type = "SE")
path <- system.file("extdata", file.path("MATS_output"), package = "maser") hypoxia <- maser(path, c("Hypoxia 0h", "Hypoxia 24h")) hypoxia_filt <- filterByCoverage(hypoxia, avg_reads = 5) hypoxia_top <- topEvents(hypoxia_filt) display(hypoxia_top, type = "SE")
Dotplot representation of splicing events.
dotplot(events, type = c("A3SS", "A5SS", "SE", "RI", "MXE"), fdr = 0.05, deltaPSI = 0.1)
dotplot(events, type = c("A3SS", "A5SS", "SE", "RI", "MXE"), fdr = 0.05, deltaPSI = 0.1)
events |
a maser object. |
type |
character indicating splice type. Possible values are
|
fdr |
numeric, FDR (False Discovery Rate) cutoff. |
deltaPSI |
numeric, absolute minimum PSI (Percent spliced-in) change |
a ggplot object.
path <- system.file("extdata", file.path("MATS_output"), package = "maser") hypoxia <- maser(path, c("Hypoxia 0h", "Hypoxia 24h")) hypoxia_filt <- filterByCoverage(hypoxia, avg_reads = 5) dotplot(hypoxia_filt, type = "SE")
path <- system.file("extdata", file.path("MATS_output"), package = "maser") hypoxia <- maser(path, c("Hypoxia 0h", "Hypoxia 24h")) hypoxia_filt <- filterByCoverage(hypoxia, avg_reads = 5) dotplot(hypoxia_filt, type = "SE")
Filter splicing events based on coverage.
filterByCoverage(events, avg_reads = 5)
filterByCoverage(events, avg_reads = 5)
events |
a maser object. |
avg_reads |
numeric, average number of reads covering the splice event. |
a maser object.
path <- system.file("extdata", file.path("MATS_output"), package = "maser") hypoxia <- maser(path, c("Hypoxia 0h", "Hypoxia 24h")) hypoxia_filt <- filterByCoverage(hypoxia, avg_reads = 5)
path <- system.file("extdata", file.path("MATS_output"), package = "maser") hypoxia <- maser(path, c("Hypoxia 0h", "Hypoxia 24h")) hypoxia_filt <- filterByCoverage(hypoxia, avg_reads = 5)
Filter splicing events based on event identifier and type.
filterByEventId(events, event_id, type = c("A3SS", "A5SS", "SE", "RI", "MXE"))
filterByEventId(events, event_id, type = c("A3SS", "A5SS", "SE", "RI", "MXE"))
events |
a maser object. |
event_id |
numeric vector of event identifiers. |
type |
character indicating splice type. Possible values are
|
a maser object.
path <- system.file("extdata", file.path("MATS_output"), package = "maser") hypoxia <- maser(path, c("Hypoxia 0h", "Hypoxia 24h")) filterByEventId(hypoxia, 33208, "SE")
path <- system.file("extdata", file.path("MATS_output"), package = "maser") hypoxia <- maser(path, c("Hypoxia 0h", "Hypoxia 24h")) filterByEventId(hypoxia, 33208, "SE")
Retrieve splicing events for a given gene.
geneEvents(events, geneS, fdr = 0.05, deltaPSI = 0.1)
geneEvents(events, geneS, fdr = 0.05, deltaPSI = 0.1)
events |
a maser object. |
geneS |
a character indicating the gene symbol. |
fdr |
numeric, FDR cutoff. |
deltaPSI |
numeric, minimum PSI change. |
a maser object.
path <- system.file("extdata", file.path("MATS_output"), package = "maser") hypoxia <- maser(path, c("Hypoxia 0h", "Hypoxia 24h")) hypoxia_mib2 <- geneEvents(hypoxia, "MIB2")
path <- system.file("extdata", file.path("MATS_output"), package = "maser") hypoxia <- maser(path, c("Hypoxia 0h", "Hypoxia 24h")) hypoxia_mib2 <- geneEvents(hypoxia, "MIB2")
Retrieve genomic ranges of splicing events from a maser object.
## S4 method for signature 'Maser' granges(x, type = c("A3SS", "A5SS", "SE", "RI", "MXE"), ...)
## S4 method for signature 'Maser' granges(x, type = c("A3SS", "A5SS", "SE", "RI", "MXE"), ...)
x |
a maser object. |
type |
a character indicating the splice type. Possible values
are |
... |
additional arguments. |
a GRangesList.
path <- system.file("extdata", file.path("MATS_output"), package = "maser") hypoxia <- maser(path, c("Hypoxia 0h", "Hypoxia 24h")) head(granges(hypoxia, type = "SE"))
path <- system.file("extdata", file.path("MATS_output"), package = "maser") hypoxia <- maser(path, c("Hypoxia 0h", "Hypoxia 24h")) head(granges(hypoxia, type = "SE"))
Mapping of splice events to UniprotKB protein features.
mapProteinFeaturesToEvents(events, tracks, by = c("feature", "category"), ncores = 1)
mapProteinFeaturesToEvents(events, tracks, by = c("feature", "category"), ncores = 1)
events |
a maser object with transcript and protein identifiers. |
tracks |
a character vector indicating valid UniprotKB features or categories. |
by |
a character vector, possible values
are |
ncores |
number of cores for multithreading (available only in OSX and Linux
machines). If Windows, |
This function performs mapping of splicing events to protein features available in the UniprotKB database. Annotation tracks of protein features mapped to the hg38 build of the human genome are retrieved from the public UniprotKB FTP. The function will overlap exons involved in the splice event with the feature genomic coordinates retrieved from UniprotKB.
Annotation can be executed either by feature or category. If categories are provided, all features within the category group will be included for annotation.
Thus, batch annotation is enabled either by using by = category
or
by providing mutilple features in the tracks
argument.
Visualization of protein features can be done
using plotUniprotKBFeatures
.
a maser object with protein feature annotation.
## Create the maser object path <- system.file("extdata", file.path("MATS_output"), package = "maser") hypoxia <- maser(path, c("Hypoxia 0h", "Hypoxia 24h")) hypoxia_filt <- filterByCoverage(hypoxia, avg_reads = 5) ## Ensembl GTF annotation for SRSF6 gtf_path <- system.file("extdata", file.path("GTF", "SRSF6_Ensembl85.gtf"), package = "maser") ens_gtf <- rtracklayer::import.gff(gtf_path) ## Retrieve gene specific splice events srsf6_events <- geneEvents(hypoxia_filt, geneS = "SRSF6") ## Map splicing events to transcripts srsf6_mapped <- mapTranscriptsToEvents(srsf6_events, ens_gtf) ## Annotate splice events with protein domains srsf6_annot <- mapProteinFeaturesToEvents(srsf6_mapped, tracks = "domain") head(annotation(srsf6_annot, "SE"))
## Create the maser object path <- system.file("extdata", file.path("MATS_output"), package = "maser") hypoxia <- maser(path, c("Hypoxia 0h", "Hypoxia 24h")) hypoxia_filt <- filterByCoverage(hypoxia, avg_reads = 5) ## Ensembl GTF annotation for SRSF6 gtf_path <- system.file("extdata", file.path("GTF", "SRSF6_Ensembl85.gtf"), package = "maser") ens_gtf <- rtracklayer::import.gff(gtf_path) ## Retrieve gene specific splice events srsf6_events <- geneEvents(hypoxia_filt, geneS = "SRSF6") ## Map splicing events to transcripts srsf6_mapped <- mapTranscriptsToEvents(srsf6_events, ens_gtf) ## Annotate splice events with protein domains srsf6_annot <- mapProteinFeaturesToEvents(srsf6_mapped, tracks = "domain") head(annotation(srsf6_annot, "SE"))
Mapping of splice events to Ensembl transcripts.
mapTranscriptsToEvents(events, gtf, ncores = 1)
mapTranscriptsToEvents(events, gtf, ncores = 1)
events |
a maser object. |
gtf |
a |
ncores |
number of cores for multithreading (available only in OSX and Linux
machines). If Windows, |
This function performs mapping of splice events in the maser object to Ensembl transcripts by overlapping exons involved in the splice event to the transcript models provided in the GTF.
Each type of splice event requires a specific mapping procedure (described below).
The mapping will also add Uniprot identifiers when the ENST transcript encodes for a protein.
Visualization of affected transcripts can be done
using plotTranscripts
.
Transcript(s) overlapping the cassette exon, as well both flanking exons (i.e upstream and downstream exons).
Transcript(s) overlapping both flanking exons but not the cassettte exon.
Transcript(s) overlapping exactly the retained intron.
Transcript(s) where intron is spliced out and overlapping both flanking exons.
Transcript(s) overlapping the first exon and both flanking exons.
Transcript(s) overlapping the second exon and both flanking exons.
Transcript(s) overlapping both short and downstream exons.
Transcript(s) overlapping both long and downstream exons.
a maser object with transcript and protein identifiers.
## Create the maser object path <- system.file("extdata", file.path("MATS_output"), package = "maser") hypoxia <- maser(path, c("Hypoxia 0h", "Hypoxia 24h")) hypoxia_filt <- filterByCoverage(hypoxia, avg_reads = 5) ## Ensembl GTF annotation for SRSF6 gtf_path <- system.file("extdata", file.path("GTF", "Ensembl85_examples.gtf.gz"), package = "maser") ens_gtf <- rtracklayer::import.gff(gtf_path) ## Retrieve gene specific splice events srsf6_events <- geneEvents(hypoxia_filt, geneS = "SRSF6") ## Map splicing events to transcripts srsf6_mapped <- mapTranscriptsToEvents(srsf6_events, ens_gtf) head(annotation(srsf6_mapped, "SE"))
## Create the maser object path <- system.file("extdata", file.path("MATS_output"), package = "maser") hypoxia <- maser(path, c("Hypoxia 0h", "Hypoxia 24h")) hypoxia_filt <- filterByCoverage(hypoxia, avg_reads = 5) ## Ensembl GTF annotation for SRSF6 gtf_path <- system.file("extdata", file.path("GTF", "Ensembl85_examples.gtf.gz"), package = "maser") ens_gtf <- rtracklayer::import.gff(gtf_path) ## Retrieve gene specific splice events srsf6_events <- geneEvents(hypoxia_filt, geneS = "SRSF6") ## Map splicing events to transcripts srsf6_mapped <- mapTranscriptsToEvents(srsf6_events, ens_gtf) head(annotation(srsf6_mapped, "SE"))
Create a maser object by importing rMATS splicing events.
maser(path, cond_labels, ftype = c("ReadsOnTargetAndJunctionCounts", "JunctionCountOnly", "JCEC", "JC"))
maser(path, cond_labels, ftype = c("ReadsOnTargetAndJunctionCounts", "JunctionCountOnly", "JCEC", "JC"))
path |
a character specifiying the folder containing rMATS output files. |
cond_labels |
a character vector of length 2 describing labels for experimental conditions. |
ftype |
a character indicating the rMATS file type.
Possible values are |
This function creates a maser object by importing rMATS output.
ftype
indicates which rMATS files to import.
ReadsOnTargetandJunction or JunctionCountOnly
are used in rMATS 3.2.5
or lower. Newer versions (>4.0.1) use "JCEC" or "JC"
nomenclature.
A maser object.
path <- system.file("extdata", file.path("MATS_output"), package = "maser") hypoxia <- maser(path, c("Hypoxia 0h", "Hypoxia 24h"))
path <- system.file("extdata", file.path("MATS_output"), package = "maser") hypoxia <- maser(path, c("Hypoxia 0h", "Hypoxia 24h"))
S4 class to represent splicing events imported from rMATS.
Prinicipal component analysis of PSI distributions.
pca(events, type = c("A3SS", "A5SS", "SE", "RI", "MXE"))
pca(events, type = c("A3SS", "A5SS", "SE", "RI", "MXE"))
events |
a maser object. |
type |
character indicating splice type. Possible values are
|
a ggplot object.
path <- system.file("extdata", file.path("MATS_output"), package = "maser") hypoxia <- maser(path, c("Hypoxia 0h", "Hypoxia 24h")) hypoxia_filt <- filterByCoverage(hypoxia, avg_reads = 5) pca(hypoxia_filt, type = "RI")
path <- system.file("extdata", file.path("MATS_output"), package = "maser") hypoxia <- maser(path, c("Hypoxia 0h", "Hypoxia 24h")) hypoxia_filt <- filterByCoverage(hypoxia, avg_reads = 5) pca(hypoxia_filt, type = "RI")
Boxplots of Percent spliced-in levels for gene events.
plotGenePSI(events, type = c("A3SS", "A5SS", "SE", "RI", "MXE"), show_replicates = TRUE)
plotGenePSI(events, type = c("A3SS", "A5SS", "SE", "RI", "MXE"), show_replicates = TRUE)
events |
a maser object. |
type |
character indicating splice type. Possible values
are |
show_replicates |
logical, add data points for individual replicates |
a ggplot object.
path <- system.file("extdata", file.path("MATS_output"), package = "maser") hypoxia <- maser(path, c("Hypoxia 0h", "Hypoxia 24h")) hypoxia_filt <- filterByCoverage(hypoxia, avg_reads = 5) hypoxia_mib2 <- geneEvents(hypoxia_filt, geneS = "MIB2") plotGenePSI(hypoxia_mib2, type = "SE", show_replicates = TRUE)
path <- system.file("extdata", file.path("MATS_output"), package = "maser") hypoxia <- maser(path, c("Hypoxia 0h", "Hypoxia 24h")) hypoxia_filt <- filterByCoverage(hypoxia, avg_reads = 5) hypoxia_mib2 <- geneEvents(hypoxia_filt, geneS = "MIB2") plotGenePSI(hypoxia_mib2, type = "SE", show_replicates = TRUE)
Mapping and visualization of Ensembl transcripts affected by splicing.
plotTranscripts(events, type = c("A3SS", "A5SS", "SE", "RI", "MXE"), event_id, gtf, zoom = FALSE, show_PSI = TRUE)
plotTranscripts(events, type = c("A3SS", "A5SS", "SE", "RI", "MXE"), event_id, gtf, zoom = FALSE, show_PSI = TRUE)
events |
a maser object. |
type |
character indicating splice type. Possible values are
|
event_id |
numeric, event identifier. |
gtf |
a |
zoom |
logical, zoom to the genomic coordinates of the splice event. |
show_PSI |
logical, display the PSI track. |
This is a wrapper function for performing both mapping and
visualization of Ensembl transcripts that are compatible with the splice
event. This function calls mapTranscriptsToEvents
for
transcript mapping, which in turn
uses findOverlaps
for transcript
overlapping. The GViz
package is used for
creating annotation tracks for genomic visualization of splicing events.
Each type of splice event requires a specific overlapping rule
(described below), #' and a customized Gviz
plot is created for
each splicing type.
Transcript(s) overlapping the cassette exon, as well both flanking exons (i.e upstream and downstream exons).
Transcript(s) overlapping both flanking exons but not the cassettte exon.
Transcript(s) overlapping exactly the retained intron.
Transcript(s) where intron is spliced out and overlapping both flanking exons.
Transcript(s) overlapping the first exon and both flanking exons.
Transcript(s) overlapping the second exon and both flanking exons.
Transcript(s) overlapping both short and downstream exons.
Transcript(s) overlapping both long and downstream exons.
a Gviz object.
## Create the maser object path <- system.file("extdata", file.path("MATS_output"), package = "maser") hypoxia <- maser(path, c("Hypoxia 0h", "Hypoxia 24h")) hypoxia_filt <- filterByCoverage(hypoxia, avg_reads = 5) ## Ensembl GTF annotation for SRSF6 gtf_path <- system.file("extdata", file.path("GTF", "SRSF6_Ensembl85.gtf"), package = "maser") ens_gtf <- rtracklayer::import.gff(gtf_path) ## Retrieve gene specific splicing events srsf6_events <- geneEvents(hypoxia_filt, geneS = "SRSF6") ## Plot exon skipping event plotTranscripts(srsf6_events, type = "SE", event_id = 33209, gtf = ens_gtf)
## Create the maser object path <- system.file("extdata", file.path("MATS_output"), package = "maser") hypoxia <- maser(path, c("Hypoxia 0h", "Hypoxia 24h")) hypoxia_filt <- filterByCoverage(hypoxia, avg_reads = 5) ## Ensembl GTF annotation for SRSF6 gtf_path <- system.file("extdata", file.path("GTF", "SRSF6_Ensembl85.gtf"), package = "maser") ens_gtf <- rtracklayer::import.gff(gtf_path) ## Retrieve gene specific splicing events srsf6_events <- geneEvents(hypoxia_filt, geneS = "SRSF6") ## Plot exon skipping event plotTranscripts(srsf6_events, type = "SE", event_id = 33209, gtf = ens_gtf)
Mapping and visualization of UniprotKB protein features affected by splicing.
plotUniprotKBFeatures(events, type = c("A3SS", "A5SS", "SE", "RI", "MXE"), event_id, gtf, features, zoom = FALSE, show_transcripts = FALSE, show_PSI = TRUE, ncores = 1)
plotUniprotKBFeatures(events, type = c("A3SS", "A5SS", "SE", "RI", "MXE"), event_id, gtf, features, zoom = FALSE, show_transcripts = FALSE, show_PSI = TRUE, ncores = 1)
events |
a maser object. |
type |
character indicating splice type. Possible values are
|
event_id |
numeric, event identifier. |
gtf |
a |
features |
a character vector indicating valid UniprotKB features. |
zoom |
logical, zoom to the genomic coordinates of the splice event. |
show_transcripts |
logical, display transcripts track. |
show_PSI |
logical, display the PSI track. |
ncores |
number of cores for multithreading (available only in OSX and Linux
machines). If Windows, |
This is a wrapper function for performing both mapping and
visualization of protein features affected by the splice event. This function
calls mapProteinFeaturesToEvents
for mapping of protein
features to splicing events.
The GViz
package is used for creating
annotation tracks for genomic visualization.
Multiple protein annotation tracks can be created using the features
argument.
a Gviz object.
## Create the maser object path <- system.file("extdata", file.path("MATS_output"), package = "maser") hypoxia <- maser(path, c("Hypoxia 0h", "Hypoxia 24h")) hypoxia_filt <- filterByCoverage(hypoxia, avg_reads = 5) ## Ensembl GTF annotation for SRSF6 gtf_path <- system.file("extdata", file.path("GTF", "SRSF6_Ensembl85.gtf"), package = "maser") ens_gtf <- rtracklayer::import.gff(gtf_path) ## Retrieve gene specific splicing events srsf6_events <- geneEvents(hypoxia_filt, geneS = "SRSF6") ## Map splicing events to transcripts srsf6_mapped <- mapTranscriptsToEvents(srsf6_events, ens_gtf) ## Plot splice event, transcripts and protein features plotUniprotKBFeatures(srsf6_mapped, "SE", event_id = 33209, gtf = ens_gtf, features = c("domain"), show_transcripts = TRUE)
## Create the maser object path <- system.file("extdata", file.path("MATS_output"), package = "maser") hypoxia <- maser(path, c("Hypoxia 0h", "Hypoxia 24h")) hypoxia_filt <- filterByCoverage(hypoxia, avg_reads = 5) ## Ensembl GTF annotation for SRSF6 gtf_path <- system.file("extdata", file.path("GTF", "SRSF6_Ensembl85.gtf"), package = "maser") ens_gtf <- rtracklayer::import.gff(gtf_path) ## Retrieve gene specific splicing events srsf6_events <- geneEvents(hypoxia_filt, geneS = "SRSF6") ## Map splicing events to transcripts srsf6_mapped <- mapTranscriptsToEvents(srsf6_events, ens_gtf) ## Plot splice event, transcripts and protein features plotUniprotKBFeatures(srsf6_mapped, "SE", event_id = 33209, gtf = ens_gtf, features = c("domain"), show_transcripts = TRUE)
Retrieve PSI (percent spliced in) values from a maser object.
PSI(events, type)
PSI(events, type)
events |
a maser object. |
type |
a character indicating the splice type. Possible values
are |
a matrix.
path <- system.file("extdata", file.path("MATS_output"), package = "maser") hypoxia <- maser(path, c("Hypoxia 0h", "Hypoxia 24h")) head(PSI(hypoxia, "SE"))
path <- system.file("extdata", file.path("MATS_output"), package = "maser") hypoxia <- maser(path, c("Hypoxia 0h", "Hypoxia 24h")) head(PSI(hypoxia, "SE"))
Retrieve PSI (percent spliced in) values from a maser object.
## S4 method for signature 'Maser,character' PSI(events, type = c("A3SS", "A5SS", "SE", "RI", "MXE"))
## S4 method for signature 'Maser,character' PSI(events, type = c("A3SS", "A5SS", "SE", "RI", "MXE"))
events |
a maser object. |
type |
a character indicating the splice type. Possible values
are |
a matrix.
path <- system.file("extdata", file.path("MATS_output"), package = "maser") hypoxia <- maser(path, c("Hypoxia 0h", "Hypoxia 24h")) head(PSI(hypoxia, "SE"))
path <- system.file("extdata", file.path("MATS_output"), package = "maser") hypoxia <- maser(path, c("Hypoxia 0h", "Hypoxia 24h")) head(PSI(hypoxia, "SE"))
Proportion of events by splicing type.
splicingDistribution(events, fdr = 0.05, deltaPSI = 0.1)
splicingDistribution(events, fdr = 0.05, deltaPSI = 0.1)
events |
a maser object. |
fdr |
numeric, FDR (False Discovery Rate) cutoff. |
deltaPSI |
numeric, absolute minimum PSI (Percent spliced-in) change |
a ggplot object.
path <- system.file("extdata", file.path("MATS_output"), package = "maser") hypoxia <- maser(path, c("Hypoxia 0h", "Hypoxia 24h")) hypoxia_filt <- filterByCoverage(hypoxia, avg_reads = 5) splicingDistribution(hypoxia_filt)
path <- system.file("extdata", file.path("MATS_output"), package = "maser") hypoxia <- maser(path, c("Hypoxia 0h", "Hypoxia 24h")) hypoxia_filt <- filterByCoverage(hypoxia, avg_reads = 5) splicingDistribution(hypoxia_filt)
Retrieve rMATS stats of differential splicing from a maser object.
## S4 method for signature 'Maser' summary(object, type = c("A3SS", "A5SS", "SE", "RI", "MXE"))
## S4 method for signature 'Maser' summary(object, type = c("A3SS", "A5SS", "SE", "RI", "MXE"))
object |
a maser object. |
type |
a character indicating the splice type. Possible values
are |
a data.frame.
path <- system.file("extdata", file.path("MATS_output"), package = "maser") hypoxia <- maser(path, c("Hypoxia 0h", "Hypoxia 24h")) head(summary(hypoxia, "SE"))
path <- system.file("extdata", file.path("MATS_output"), package = "maser") hypoxia <- maser(path, c("Hypoxia 0h", "Hypoxia 24h")) head(summary(hypoxia, "SE"))
Filter splicing events based on false discovery rate and PSI change.
topEvents(events, fdr = 0.05, deltaPSI = 0.1)
topEvents(events, fdr = 0.05, deltaPSI = 0.1)
events |
a maser object. |
fdr |
numeric, FDR (False Discovery Rate) cutoff. |
deltaPSI |
numeric, absolute minimum PSI (Percent spliced-in) change |
a maser object.
path <- system.file("extdata", file.path("MATS_output"), package = "maser") hypoxia <- maser(path, c("Hypoxia 0h", "Hypoxia 24h")) ## To select all events with minimum 10% change in PSI, and FDR < 0.01 hypoxia_top <- topEvents(hypoxia, fdr = 0.01, deltaPSI = 0.1)
path <- system.file("extdata", file.path("MATS_output"), package = "maser") hypoxia <- maser(path, c("Hypoxia 0h", "Hypoxia 24h")) ## To select all events with minimum 10% change in PSI, and FDR < 0.01 hypoxia_top <- topEvents(hypoxia, fdr = 0.01, deltaPSI = 0.1)
Volcano plot of splicing events.
volcano(events, type = c("A3SS", "A5SS", "SE", "RI", "MXE"), fdr = 0.05, deltaPSI = 0.1)
volcano(events, type = c("A3SS", "A5SS", "SE", "RI", "MXE"), fdr = 0.05, deltaPSI = 0.1)
events |
a maser object. |
type |
character indicating splice type. Possible values are
|
fdr |
numeric, FDR (False Discovery Rate) cutoff. |
deltaPSI |
numeric, absolute minimum PSI (Percent spliced-in) change |
a ggplot object.
path <- system.file("extdata", file.path("MATS_output"), package = "maser") hypoxia <- maser(path, c("Hypoxia 0h", "Hypoxia 24h")) hypoxia_filt <- filterByCoverage(hypoxia, avg_reads = 5) volcano(hypoxia_filt, type = "SE")
path <- system.file("extdata", file.path("MATS_output"), package = "maser") hypoxia <- maser(path, c("Hypoxia 0h", "Hypoxia 24h")) hypoxia_filt <- filterByCoverage(hypoxia, avg_reads = 5) volcano(hypoxia_filt, type = "SE")