Move one or more columns from a GRangesList elements into assays in a RangesSummarizedEperiment

grlToSE(x, ...)

# S4 method for class 'GRangesList'
grlToSE(
  x,
  assayCols = c(),
  metaCols = c(),
  keyvals = c(),
  by = c("min", "max"),
  ...,
  ignore.strand = FALSE
)

Arguments

x

A GrangesList

...

Passed to reduce

assayCols

Columns to move to separate assays

metaCols

Columns to move to mcols within the rowRanges element

keyvals

The value to use when choosing representative values

by

How to choose by keyvals

ignore.strand

logical(1). Whether the strand of the input ranges should be ignored or not.

Value

A RangedSummarizedExperiment

Details

Given a GRangesList which would commonly represent multiple samples, reduce any overlapping ranges into a consensus range, setting any metadata columns to be retained as separate assays. These columns may contain values such as coverage, p-values etc.

Additional columns can also be placed as rowData columns where the original values are better suited to information about the consensus range rather than the sample (or GRangesList element).

Only one value for each range will be retained, and these are chosen using the value provided as the keyvals, taking either the min or max value in this column as the representative range.

Examples

a <- GRanges("chr1:1-10")
a$feature <- "Gene"
a$p <- 0.1
b <- GRanges(c("chr1:6-15", "chr1:15"))
b$feature <- c("Gene", "Promoter")
b$p <- c(0.5, 0.01)
grl <- GRangesList(a = a, b = b)
grl
#> GRangesList object of length 2:
#> $a
#> GRanges object with 1 range and 2 metadata columns:
#>       seqnames    ranges strand |     feature         p
#>          <Rle> <IRanges>  <Rle> | <character> <numeric>
#>   [1]     chr1      1-10      * |        Gene       0.1
#>   -------
#>   seqinfo: 1 sequence from an unspecified genome; no seqlengths
#> 
#> $b
#> GRanges object with 2 ranges and 2 metadata columns:
#>       seqnames    ranges strand |     feature         p
#>          <Rle> <IRanges>  <Rle> | <character> <numeric>
#>   [1]     chr1      6-15      * |        Gene      0.50
#>   [2]     chr1        15      * |    Promoter      0.01
#>   -------
#>   seqinfo: 1 sequence from an unspecified genome; no seqlengths
#> 
se <- grlToSE(
  grl, assayCols = "p", metaCols = "feature", keyvals = "p", by = "min"
)
assay(se, "p")
#>        a    b
#> [1,] 0.1 0.01
rowRanges(se)
#> GRanges object with 1 range and 1 metadata column:
#>       seqnames    ranges strand |         feature
#>          <Rle> <IRanges>  <Rle> | <CharacterList>
#>   [1]     chr1      1-15      * |   Gene,Promoter
#>   -------
#>   seqinfo: 1 sequence from an unspecified genome; no seqlengths