Last data update: 2014.03.03

R: BEDFile objects
BEDFile-classR Documentation

BEDFile objects

Description

These functions support the import and export of the UCSC BED format and its variants, including BEDGraph.

Usage

## S4 method for signature 'BEDFile,ANY,ANY'
import(con, format, text, trackLine = TRUE,
                   genome = NA, colnames = NULL,
                   which = NULL, seqinfo = NULL, extraCols = character())
import.bed(con, ...)
import.bed15(con, ...)
import.bedGraph(con,  ...)

## S4 method for signature 'ANY,BEDFile,ANY'
export(object, con, format, ...)
## S4 method for signature 'GenomicRanges,BEDFile,ANY'
export(object, con, format,
                  append = FALSE, index = FALSE,
                  ignore.strand = FALSE, trackLine = NULL)
## S4 method for signature 'UCSCData,BEDFile,ANY'
export(object, con, format,
                   trackLine = TRUE, ...)
export.bed(object, con, ...)
export.bed15(object, con, ...)
## S4 method for signature 'GenomicRanges,BED15File,ANY'
export(object, con, format,
                  expNames = NULL, trackLine = NULL, ...)
export.bedGraph(object, con, ...)

Arguments

con

A path, URL, connection or BEDFile object. For the functions ending in .bed, .bedGraph and .bed15, the file format is indicated by the function name. For the base export and import functions, the format must be indicated another way. If con is a path, URL or connection, either the file extension or the format argument needs to be one of “bed”, “bed15”, “bedGraph” or “bedpe”. Compressed files (“gz”, “bz2” and “xz”) are handled transparently.

object

The object to export, should be a GRanges or something coercible to a GRanges. If the object has a method for asBED (like GRangesList), it is called prior to coercion. This makes it possible to export a GRangesList or TxDb in a way that preserves the hierarchical structure. For exporting multiple tracks, in the UCSC track line metaformat, pass a GenomicRangesList, or something coercible to one.

format

If not missing, should be one of “bed”, “bed15”, “bedGraph” or “bedpe”.

text

If con is missing, a character vector to use as the input

trackLine

For import, an imported track line will be stored in a TrackLine object, as part of the returned UCSCData. For the UCSCData method on export, whether to output the UCSC track line stored on the object, for the other export methods, the actual TrackLine object to export.

genome

The identifier of a genome, or NA if unknown. Typically, this is a UCSC identifier like “hg19”. An attempt will be made to derive the seqinfo on the return value using either an installed BSgenome package or UCSC, if network access is available.

colnames

A character vector naming the columns to parse. These should name columns in the result, not those in the BED spec, so e.g. specify “thick”, instead of “thickStart”.

which

A GRanges or other range-based object supported by findOverlaps. Only the intervals in the file overlapping the given ranges are returned. This is much more efficient when the file is indexed with the tabix utility.

index

If TRUE, automatically compress and index the output file with bgzf and tabix. Note that tabix indexing will sort the data by chromosome and start. Tabix supports a single track in a file.

ignore.strand

Whether to output the strand when not required (by the existence of later fields).

seqinfo

If not NULL, the Seqinfo object to set on the result. If the genome argument is not NA, it must agree with genome(seqinfo).

extraCols

A character vector in the same form as colClasses from read.table. It should indicate the name and class of each extra/special column to read from the BED file. As BED does not encode column names, these are assumed to be the last columns in the file. This enables parsing of the various BEDX+Y formats.

append

If TRUE, and con points to a file path, the data is appended to the file. Obviously, if con is a connection, the data is always appended.

expNames

character vector of column names in object to export as sample columns in the BED15 file.

...

Arguments to pass down to methods to other methods. For import, the flow eventually reaches the BEDFile method on import. When trackLine is TRUE or the target format is BED15, the arguments are passed through export.ucsc, so track line parameters are supported.

Details

The BED format is a tab-separated table of intervals, with annotations like name, score and even sub-intervals for representing alignments and gene models. Official (UCSC) child formats currently include BED15 (adding a number matrix for e.g. expression data across multiple samples) and BEDGraph (a compressed means of storing a single score variable, e.g. coverage; overlapping features are not allowed). Many tools and organizations have extended the BED format with additional columns for particular use cases. Currently only BEDPE is supported by rtracklayer; others are not yet supported but a mechanism will be added soon. The advantage of BED is its balance between simplicity and expressiveness. It is also relatively scalable, because only the first three columns (chrom, start and end) are required. Thus, BED is best suited for representing simple features. For specialized cases, one is usually better off with another format. For example, genome-scale vectors belong in BigWig, alignments from high-throughput sequencing belong in BAM, and gene models are more richly expressed in GFF.

The following is the mapping of BED elements to a GRanges object. NA values are allowed only where indicated. These appear as a “.” in the file. Only the first three columns (chrom, start and strand) are required. The other columns can only be included if all previous columns (to the left) are included. Upon export, default values are used to automatically pad the table, if necessary.

chrom, start, end

the ranges component.

name

character vector (NA's allowed) in the name column; defaults to NA on export.

score

numeric vector in the score column, accessible via the score accessor. Defaults to 0 on export. This is the only column present in BEDGraph (besides chrom, start and end), and it is required.

strand

strand factor (NA's allowed) in the strand column, accessible via the strand accessor; defaults to NA on export.

thickStart, thickEnd

Ranges object in a column named thick; defaults to the ranges of the feature on export.

itemRgb

an integer matrix of color codes, as returned by col2rgb, or any valid input to col2rgb, in the itemRgb column; default is NA on export, which translates to black.

blockSizes, blockStarts, blockCounts

RangesList object in a column named blocks; defaults to empty upon BED15 export.

These columns are present only in BED15:

expCount, expIds, expScores

A column for each unique element in expIds, containing the corresponding values from expScores. When a value is not present for a feature, NA is substituted. NA values become -10000 in the file.

Value

For a “bedpe” file, a Pairs object combining two GRanges. The name and score are carried over to the metadata columns.

Otherwise, a GRanges with the metadata columns described in the details.

BEDX+Y formats

To import one of the multitude of BEDX+Y formats, such as those used to distribute ENCODE data through UCSC (narrowPeaks, etc), specify the extraCols argument to indicate the expected names and classes of the special columns. We assume that the last length(extraCols) columns are special, and that the preceding columns adhere to the BED format.

BEDFile objects

The BEDFile class extends RTLFile and is a formal represention of a resource in the BED format. To cast a path, URL or connection to a BEDFile, pass it to the BEDFile constructor. Classes and constructors also exist for the subclasses BED15File, BEDGraphFile and BEDPEFile.

Author(s)

Michael Lawrence

References

http://genome.ucsc.edu/goldenPath/help/customTrack.html http://bedtools.readthedocs.org/en/latest/content/general-usage.html

Examples

  test_path <- system.file("tests", package = "rtracklayer")
  test_bed <- file.path(test_path, "test.bed")

  test <- import(test_bed)
  test
  
  test_bed_file <- BEDFile(test_bed)
  import(test_bed_file)
  
  test_bed_con <- file(test_bed)
  import(test_bed_con, format = "bed")
  close(test_bed_con)
  
  import(test_bed, trackLine = FALSE)
  import(test_bed, genome = "hg19")
  import(test_bed, colnames = c("name", "strand", "thick"))

  which <- GRanges("chr7:1-127473000")
  import(test_bed, which = which)

## Not run: 
  test_bed_out <- file.path(tempdir(), "test.bed")
  export(test, test_bed_out)
  
  test_bed_out_file <- BEDFile(test_bed_out)
  export(test, test_bed_out_file)

  export(test, test_bed_out, name = "Alternative name")
  
  test_bed_gz <- paste(test_bed_out, ".gz", sep = "")
  export(test, test_bed_gz)
  
  export(test, test_bed_out, index = TRUE)
  export(test, test_bed_out, index = TRUE, trackLine = FALSE)
  
  bed_text <- export(test, format = "bed")
  test <- import(format = "bed", text = bed_text)

## End(Not run)

Results


R version 3.3.1 (2016-06-21) -- "Bug in Your Hair"
Copyright (C) 2016 The R Foundation for Statistical Computing
Platform: x86_64-pc-linux-gnu (64-bit)

R is free software and comes with ABSOLUTELY NO WARRANTY.
You are welcome to redistribute it under certain conditions.
Type 'license()' or 'licence()' for distribution details.

R is a collaborative project with many contributors.
Type 'contributors()' for more information and
'citation()' on how to cite R or R packages in publications.

Type 'demo()' for some demos, 'help()' for on-line help, or
'help.start()' for an HTML browser interface to help.
Type 'q()' to quit R.

> library(rtracklayer)
Loading required package: GenomicRanges
Loading required package: BiocGenerics
Loading required package: parallel

Attaching package: 'BiocGenerics'

The following objects are masked from 'package:parallel':

    clusterApply, clusterApplyLB, clusterCall, clusterEvalQ,
    clusterExport, clusterMap, parApply, parCapply, parLapply,
    parLapplyLB, parRapply, parSapply, parSapplyLB

The following objects are masked from 'package:stats':

    IQR, mad, xtabs

The following objects are masked from 'package:base':

    Filter, Find, Map, Position, Reduce, anyDuplicated, append,
    as.data.frame, cbind, colnames, do.call, duplicated, eval, evalq,
    get, grep, grepl, intersect, is.unsorted, lapply, lengths, mapply,
    match, mget, order, paste, pmax, pmax.int, pmin, pmin.int, rank,
    rbind, rownames, sapply, setdiff, sort, table, tapply, union,
    unique, unsplit

Loading required package: S4Vectors
Loading required package: stats4

Attaching package: 'S4Vectors'

The following objects are masked from 'package:base':

    colMeans, colSums, expand.grid, rowMeans, rowSums

Loading required package: IRanges
Loading required package: GenomeInfoDb
> png(filename="/home/ddbj/snapshot/RGM3/R_BC/result/rtracklayer/BEDFile-class.Rd_%03d_medium.png", width=480, height=480)
> ### Name: BEDFile-class
> ### Title: BEDFile objects
> ### Aliases: class:BEDFile class:BED15File class:BEDGraphFile BEDFile-class
> ###   BED15File-class BEDGraphFile-class BEDPEFile-class BEDFile BED15File
> ###   BEDGraphFile BEDPEFile import,BEDFile,ANY,ANY-method
> ###   import,BED15File,ANY,ANY-method import,BEDPEFile,ANY,ANY-method
> ###   import.bed import.bed,ANY-method import.bed15 import.bed15,ANY-method
> ###   import.bedGraph import.bedGraph,ANY-method
> ###   export,ANY,BEDFile,ANY-method export,GenomicRanges,BEDFile,ANY-method
> ###   export,GRangesList,BEDFile,ANY-method
> ###   export,GenomicRangesList,BEDFile,ANY-method
> ###   export,UCSCData,BEDFile,ANY-method
> ###   export,GenomicRanges,BED15File,ANY-method export.bed
> ###   export.bed,ANY-method export.bed15 export.bed15,ANY-method
> ###   export.bedGraph export.bedGraph,ANY-method
> ### Keywords: methods classes
> 
> ### ** Examples
> 
>   test_path <- system.file("tests", package = "rtracklayer")
>   test_bed <- file.path(test_path, "test.bed")
> 
>   test <- import(test_bed)
>   test
UCSC track 'ItemRGBDemo'
UCSCData object with 5 ranges and 5 metadata columns:
      seqnames                 ranges strand |        name     score
         <Rle>              <IRanges>  <Rle> | <character> <numeric>
  [1]     chr7 [127471197, 127472363]      + |        Pos1         0
  [2]     chr7 [127472364, 127473530]      + |        Pos2         2
  [3]     chr7 [127473531, 127474697]      - |        Neg1         0
  [4]     chr9 [127474698, 127475864]      + |        Pos3         5
  [5]     chr9 [127475865, 127477031]      - |        Neg2         5
          itemRgb                  thick                                 blocks
      <character>              <IRanges>                          <IRangesList>
  [1]     #FF0000 [127471197, 127472363] [   1,  300] [ 501,  700] [1068, 1167]
  [2]     #FF0000 [127472364, 127473530]                [  1,  250] [668, 1167]
  [3]     #FF0000 [127473531, 127474697]                              [1, 1167]
  [4]     #FF0000 [127474698, 127475864]                              [1, 1167]
  [5]     #0000FF [127475865, 127477031]                              [1, 1167]
  -------
  seqinfo: 2 sequences from an unspecified genome; no seqlengths
>   
>   test_bed_file <- BEDFile(test_bed)
>   import(test_bed_file)
UCSC track 'ItemRGBDemo'
UCSCData object with 5 ranges and 5 metadata columns:
      seqnames                 ranges strand |        name     score
         <Rle>              <IRanges>  <Rle> | <character> <numeric>
  [1]     chr7 [127471197, 127472363]      + |        Pos1         0
  [2]     chr7 [127472364, 127473530]      + |        Pos2         2
  [3]     chr7 [127473531, 127474697]      - |        Neg1         0
  [4]     chr9 [127474698, 127475864]      + |        Pos3         5
  [5]     chr9 [127475865, 127477031]      - |        Neg2         5
          itemRgb                  thick                                 blocks
      <character>              <IRanges>                          <IRangesList>
  [1]     #FF0000 [127471197, 127472363] [   1,  300] [ 501,  700] [1068, 1167]
  [2]     #FF0000 [127472364, 127473530]                [  1,  250] [668, 1167]
  [3]     #FF0000 [127473531, 127474697]                              [1, 1167]
  [4]     #FF0000 [127474698, 127475864]                              [1, 1167]
  [5]     #0000FF [127475865, 127477031]                              [1, 1167]
  -------
  seqinfo: 2 sequences from an unspecified genome; no seqlengths
>   
>   test_bed_con <- file(test_bed)
>   import(test_bed_con, format = "bed")
UCSC track 'ItemRGBDemo'
UCSCData object with 5 ranges and 5 metadata columns:
      seqnames                 ranges strand |        name     score
         <Rle>              <IRanges>  <Rle> | <character> <numeric>
  [1]     chr7 [127471197, 127472363]      + |        Pos1         0
  [2]     chr7 [127472364, 127473530]      + |        Pos2         2
  [3]     chr7 [127473531, 127474697]      - |        Neg1         0
  [4]     chr9 [127474698, 127475864]      + |        Pos3         5
  [5]     chr9 [127475865, 127477031]      - |        Neg2         5
          itemRgb                  thick                                 blocks
      <character>              <IRanges>                          <IRangesList>
  [1]     #FF0000 [127471197, 127472363] [   1,  300] [ 501,  700] [1068, 1167]
  [2]     #FF0000 [127472364, 127473530]                [  1,  250] [668, 1167]
  [3]     #FF0000 [127473531, 127474697]                              [1, 1167]
  [4]     #FF0000 [127474698, 127475864]                              [1, 1167]
  [5]     #0000FF [127475865, 127477031]                              [1, 1167]
  -------
  seqinfo: 2 sequences from an unspecified genome; no seqlengths
>   close(test_bed_con)
>   
>   import(test_bed, trackLine = FALSE)
GRanges object with 5 ranges and 5 metadata columns:
      seqnames                 ranges strand |        name     score
         <Rle>              <IRanges>  <Rle> | <character> <numeric>
  [1]     chr7 [127471197, 127472363]      + |        Pos1         0
  [2]     chr7 [127472364, 127473530]      + |        Pos2         2
  [3]     chr7 [127473531, 127474697]      - |        Neg1         0
  [4]     chr9 [127474698, 127475864]      + |        Pos3         5
  [5]     chr9 [127475865, 127477031]      - |        Neg2         5
          itemRgb                  thick                                 blocks
      <character>              <IRanges>                          <IRangesList>
  [1]     #FF0000 [127471197, 127472363] [   1,  300] [ 501,  700] [1068, 1167]
  [2]     #FF0000 [127472364, 127473530]                [  1,  250] [668, 1167]
  [3]     #FF0000 [127473531, 127474697]                              [1, 1167]
  [4]     #FF0000 [127474698, 127475864]                              [1, 1167]
  [5]     #0000FF [127475865, 127477031]                              [1, 1167]
  -------
  seqinfo: 2 sequences from an unspecified genome; no seqlengths
>   import(test_bed, genome = "hg19")
UCSC track 'ItemRGBDemo'
UCSCData object with 5 ranges and 5 metadata columns:
      seqnames                 ranges strand |        name     score
         <Rle>              <IRanges>  <Rle> | <character> <numeric>
  [1]     chr7 [127471197, 127472363]      + |        Pos1         0
  [2]     chr7 [127472364, 127473530]      + |        Pos2         2
  [3]     chr7 [127473531, 127474697]      - |        Neg1         0
  [4]     chr9 [127474698, 127475864]      + |        Pos3         5
  [5]     chr9 [127475865, 127477031]      - |        Neg2         5
          itemRgb                  thick                                 blocks
      <character>              <IRanges>                          <IRangesList>
  [1]     #FF0000 [127471197, 127472363] [   1,  300] [ 501,  700] [1068, 1167]
  [2]     #FF0000 [127472364, 127473530]                [  1,  250] [668, 1167]
  [3]     #FF0000 [127473531, 127474697]                              [1, 1167]
  [4]     #FF0000 [127474698, 127475864]                              [1, 1167]
  [5]     #0000FF [127475865, 127477031]                              [1, 1167]
  -------
  seqinfo: 93 sequences (1 circular) from hg19 genome
>   import(test_bed, colnames = c("name", "strand", "thick"))
UCSC track 'ItemRGBDemo'
UCSCData object with 5 ranges and 2 metadata columns:
      seqnames                 ranges strand |        name
         <Rle>              <IRanges>  <Rle> | <character>
  [1]     chr7 [127471197, 127472363]      + |        Pos1
  [2]     chr7 [127472364, 127473530]      + |        Pos2
  [3]     chr7 [127473531, 127474697]      - |        Neg1
  [4]     chr9 [127474698, 127475864]      + |        Pos3
  [5]     chr9 [127475865, 127477031]      - |        Neg2
                       thick
                   <IRanges>
  [1] [127471197, 127472363]
  [2] [127472364, 127473530]
  [3] [127473531, 127474697]
  [4] [127474698, 127475864]
  [5] [127475865, 127477031]
  -------
  seqinfo: 2 sequences from an unspecified genome; no seqlengths
> 
>   which <- GRanges("chr7:1-127473000")
>   import(test_bed, which = which)
UCSC track 'ItemRGBDemo'
UCSCData object with 2 ranges and 5 metadata columns:
      seqnames                 ranges strand |        name     score
         <Rle>              <IRanges>  <Rle> | <character> <numeric>
  [1]     chr7 [127471197, 127472363]      + |        Pos1         0
  [2]     chr7 [127472364, 127473530]      + |        Pos2         2
          itemRgb                  thick                                 blocks
      <character>              <IRanges>                          <IRangesList>
  [1]     #FF0000 [127471197, 127472363] [   1,  300] [ 501,  700] [1068, 1167]
  [2]     #FF0000 [127472364, 127473530]                [  1,  250] [668, 1167]
  -------
  seqinfo: 2 sequences from an unspecified genome; no seqlengths
> 
> ## Not run: 
> ##D   test_bed_out <- file.path(tempdir(), "test.bed")
> ##D   export(test, test_bed_out)
> ##D   
> ##D   test_bed_out_file <- BEDFile(test_bed_out)
> ##D   export(test, test_bed_out_file)
> ##D 
> ##D   export(test, test_bed_out, name = "Alternative name")
> ##D   
> ##D   test_bed_gz <- paste(test_bed_out, ".gz", sep = "")
> ##D   export(test, test_bed_gz)
> ##D   
> ##D   export(test, test_bed_out, index = TRUE)
> ##D   export(test, test_bed_out, index = TRUE, trackLine = FALSE)
> ##D   
> ##D   bed_text <- export(test, format = "bed")
> ##D   test <- import(format = "bed", text = bed_text)
> ## End(Not run)
> 
> 
> 
> 
> 
> dev.off()
null device 
          1 
>