BiomartGeneRegionTrack-class | R Documentation |
A class to hold gene model data for a genomic region fetched dynamically from EBI's Biomart Ensembl data source.
## S4 method for signature 'BiomartGeneRegionTrack'
initialize(
.Object,
start = NULL,
end = NULL,
biomart,
filter = list(),
range,
genome = NULL,
chromosome = NULL,
strand = NULL,
featureMap = NULL,
symbol = NULL,
gene = NULL,
transcript = NULL,
entrez = NULL,
...
)
BiomartGeneRegionTrack(
start = NULL,
end = NULL,
biomart,
chromosome = NULL,
strand,
genome = NULL,
stacking = "squish",
filters = list(),
featureMap = NULL,
name = "BiomartGeneRegionTrack",
symbol = NULL,
gene = NULL,
entrez = NULL,
transcript = NULL,
...
)
## S4 method for signature 'BiomartGeneRegionTrack'
subset(x, from, to, chromosome, use.defaults = TRUE, ...)
.Object |
.Object |
start |
An integer scalar with the genomic start coordinates for the gene model range. |
end |
An integer scalar with the genomic end coordinates for the gene model range. |
biomart |
An optional |
filter |
filter |
range |
range |
genome |
The genome on which the track's ranges are defined. Usually
this is a valid UCSC genome identifier, however this is not being formally
checked at this point. If no mapping from genome to Biomart Ensembl data
source is possible, the |
chromosome |
The chromosome on which the track's genomic ranges are
defined. A valid UCSC chromosome identifier. Please note that at this stage
only syntactic checking takes place, i.e., the argument value needs to be a
single integer, numeric character or a character of the form |
strand |
Character scalar, the strand for which to fetch gene
information from Biomart. One in |
featureMap |
Named character vector or list to map between the fields in the Biomart data base and the features as they are used to construct the track. If multiple values are provided in a single list item, the package will use the first one that is defined in the selected Biomart. |
symbol, transcript, gene, entrez |
Character vector giving one or several gene symbols, Ensembl transcript identifiers, Ensembl gene identifiers, or ENTREZ gene identifiers, respectively. The genomic locus of their gene model will be fetch from Biomart instead of providing explicit start and end coordinates. |
... |
Additional items which will all be interpreted as further
display parameters. See |
stacking |
The stacking type for overlapping items of the track. One in
|
filters |
A list of additional filters to be applied in the Biomart
query. See |
name |
Character scalar of the track's name used in the title panel when plotting. |
x |
A valid track object class name, or the object itself, in which case the class is derived directly from it. |
from, to |
from, to |
use.defaults |
|
A track containing all gene models in a particular region as fetched from
EBI's Biomart service. Usually the user does not have to take care of the
Biomart connection, which will be established automatically based on the
provided genome and chromosome information. However, for full flexibility a
valid Mart
object may be passed on to the constructor.
Please note that this assumes a connection to one of the Ensembl gene data
sources, mapping the available query data back to the internal object slots.
The return value of the constructor function is a new object of class
BiomartGeneRegionTrack
.
initialize(BiomartGeneRegionTrack)
: Initialize.
BiomartGeneRegionTrack()
: Constructor function for
BiomartGeneRegionTrack-class
.
subset(BiomartGeneRegionTrack)
: subset a BiomartGeneRegionTrack
by coordinates and sort if necessary.
Objects can be created using the constructor function
BiomartGeneRegionTrack
.
Florian Hahne
EBI Biomart webservice at http://www.biomart.org.
DisplayPars
GdObject
GRanges
HighlightTrack
ImageMap
IRanges
RangeTrack
DataTrack
collapsing
grouping
panel.grid
plotTracks
settings
## Construct the object
## Not run:
bmTrack <- BiomartGeneRegionTrack(
start = 26682683, end = 26711643,
chromosome = 7, genome = "mm9"
)
## End(Not run)
## Plotting
plotTracks(bmTrack)
## Track names
names(bmTrack)
names(bmTrack) <- "foo"
plotTracks(bmTrack)
## Subsetting and splitting
subTrack <- subset(bmTrack, from = 26700000, to = 26705000)
length(subTrack)
subTrack <- bmTrack[transcript(bmTrack) == "ENSMUST00000144140"]
split(bmTrack, transcript(bmTrack))
## Accessors
start(bmTrack)
end(bmTrack)
width(bmTrack)
position(bmTrack)
width(subTrack) <- width(subTrack) + 100
strand(bmTrack)
strand(subTrack) <- "-"
chromosome(bmTrack)
chromosome(subTrack) <- "chrX"
genome(bmTrack)
genome(subTrack) <- "hg19"
range(bmTrack)
ranges(bmTrack)
## Annotation
identifier(bmTrack)
identifier(bmTrack, "lowest")
identifier(subTrack) <- "bar"
feature(bmTrack)
feature(subTrack) <- "foo"
exon(bmTrack)
exon(subTrack) <- letters[1:2]
gene(bmTrack)
gene(subTrack) <- "bar"
symbol(bmTrack)
symbol(subTrack) <- "foo"
transcript(bmTrack)
transcript(subTrack) <- c("foo", "bar")
chromosome(subTrack) <- "chr7"
plotTracks(subTrack)
values(bmTrack)
## Grouping
group(bmTrack)
group(subTrack) <- "Group 1"
transcript(subTrack)
plotTracks(subTrack)
## Stacking
stacking(bmTrack)
stacking(bmTrack) <- "dense"
plotTracks(bmTrack)
## coercion
as(bmTrack, "data.frame")
as(bmTrack, "UCSCData")
## HTML image map
coords(bmTrack)
tags(bmTrack)
bmTrack <- plotTracks(bmTrack)$foo
coords(bmTrack)
tags(bmTrack)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.