import.gff {rtracklayer} | R Documentation |
These are the functions for importing
RangedData
instances from connections or text.
import.gff(con, version = c("1", "2", "3"), genome = NULL, asRangedData = TRUE, colnames = NULL) import.gff1(con, ...) import.gff2(con, ...) import.gff3(con, ...) import.bed(con, variant = c("base", "bedGraph", "bed15"), trackLine = TRUE, genome = NULL, asRangedData = TRUE, colnames = NULL, ...) import.bed15(con, genome = NULL, asRangedData = TRUE, ...) import.bedGraph(con, genome = NULL, asRangedData = TRUE, ...) import.wig(con, genome = NULL, asRangedData = TRUE, ...) import.ucsc(con, subformat = c("auto", "gff1", "wig", "bed", "bed15", "bedGraph"), drop = FALSE, asRangedData = TRUE, ...) ## not yet supported on Windows import.bw(con, ...)
con |
The connection, filename or URL from which to receive the input. |
version |
The version of GFF ("1", "2" or "3"). |
genome |
The genome to set on the imported track. |
asRangedData |
A logical value. If |
variant |
Variant of BED lines, not for the user. |
trackLine |
Whether the BED data has a track line (it normally does though track lines are not mandatory). |
subformat |
The expected subformat of the UCSC data. If "auto", automatic detection of the subformat is attempted. |
drop |
If |
colnames |
Character vector indicating which columns (excluding
the required sequence name, start and end) should be imported. If
|
... |
For |
For all but import.ucsc
, an instance of
RangedData
(or one of its
subclasses) or GRanges
if asRangedData
is TRUE
or FALSE
respectively.
For import.ucsc
when drop
is FALSE
, an instance
of RangedDataList
or
GRangesList
if asRangedData
is TRUE
or FALSE
respectively.
Michael Lawrence and Patrick Aboyoun
http://genome.ucsc.edu/goldenPath/help/customTrack.html\#BED
import
for the high-level interface to these routines.
# import a GFF V2 file gffRD <- import.gff(system.file("tests", "v2.gff", package = "rtracklayer"), version = "2") gffGR <- import.gff(system.file("tests", "v2.gff", package = "rtracklayer"), version = "2", asRangedData = FALSE) # or gffRD <- import.gff2(system.file("tests", "v2.gff", package = "rtracklayer")) gffGR <- import.gff2(system.file("tests", "v2.gff", package = "rtracklayer"), asRangedData = FALSE) # import a WIG file wigRD <- import.wig(system.file("tests", "bed.wig", package = "rtracklayer")) wigGR <- import.wig(system.file("tests", "bed.wig", package = "rtracklayer"), asRangedData = FALSE) # or wigRD <- import.ucsc(system.file("tests", "bed.wig", package = "rtracklayer"), subformat = "wig", drop = TRUE) wigGR <- import.ucsc(system.file("tests", "bed.wig", package = "rtracklayer"), subformat = "wig", drop = TRUE, asRangedData = FALSE) # bigWig ## Not run: bw <- import(system.file("tests", "test.bw", package = "rtracklayer"), ranges = GenomicRanges::GRanges("chr19", IRanges(1, 6e7))) ## End(Not run)