Convert multiple Genomic objects to tibbles
# S3 method for class 'DataFrame'
as_tibble(x, rangeAsChar = TRUE, ...)
# S3 method for class 'GenomicRanges'
as_tibble(x, rangeAsChar = TRUE, name = "range", ...)
# S3 method for class 'Seqinfo'
as_tibble(x, ...)
# S3 method for class 'GInteractions'
as_tibble(x, rangeAsChar = TRUE, suffix = c(".x", ".y"), ...)
# S3 method for class 'SummarizedExperiment'
as_tibble(x, rangeAsChar = TRUE, ...)
# S3 method for class 'TopTags'
as_tibble(x, ...)
A Genomic Ranges or DataFrame object
Convert any GRanges element to a character vector
Passed to tibble::as_tibble()
Name of column to use for ranges. Ignored if rangeAsChar =
FALSE
Suffix appended to column names for anchor1 and anchor2 of a GInteractions object. Only used if specifying rangeAsChar = FALSE
A tibble
Quick and dirty conversion into a tibble.
By default, GenomicRanges will be returned with the range as a character
column called range
and all mcols parsed as the remaining columns.
Seqinfo information will be lost during coercion.
Given that names for ranges are considered as rownames in the mcols element,
these can be simply parsed by setting rownames = "id"
in the call to
as_tibble()
When coercing a DataFrame, any Compressed/SimpleList columns will be coerced to S3 list columns. Any GRanges columns will be returned as a character column, losing any additional mcols from these secondary ranges
Coercion of SummarizedExperiment objects will be performed on the
rowRanges()
element, whilst for a GInteractions
object, both anchors will returned with the default suffixes .x
and .y
Defined as an S3 method for consistency with existing tidy methods
gr <- GRanges("chr1:1-10")
gr$p_value <- runif(1)
names(gr) <- "range1"
gr
#> GRanges object with 1 range and 1 metadata column:
#> seqnames ranges strand | p_value
#> <Rle> <IRanges> <Rle> | <numeric>
#> range1 chr1 1-10 * | 0.718335
#> -------
#> seqinfo: 1 sequence from an unspecified genome; no seqlengths
as_tibble(gr)
#> # A tibble: 1 × 2
#> range p_value
#> <chr> <dbl>
#> 1 chr1:1-10 0.718
as_tibble(gr, rownames = "id")
#> # A tibble: 1 × 3
#> range id p_value
#> <chr> <chr> <dbl>
#> 1 chr1:1-10 range1 0.718
as_tibble(mcols(gr))
#> # A tibble: 1 × 1
#> p_value
#> <dbl>
#> 1 0.718
as_tibble(seqinfo(gr))
#> # A tibble: 1 × 4
#> seqnames seqlengths is_circular genome
#> <chr> <int> <lgl> <chr>
#> 1 chr1 NA NA NA
hic <- InteractionSet::GInteractions(gr, GRanges("chr1:201-210"))
hic$id <- "interaction1"
as_tibble(hic)
#> # A tibble: 1 × 4
#> anchor1 anchor2 anchor1.p_value id
#> <chr> <chr> <dbl> <chr>
#> 1 chr1:1-10 chr1:201-210 0.718 interaction1