R/mapGrlCols.R
mapGrlCols.Rd
Make consensus peaks and add individual columns from each original GRangesList element
mapGrlCols(
x,
var = NULL,
collapse = NULL,
collapse_sep = "; ",
name_sep = "_",
include = FALSE,
...
)
GRangesList
Column(s) to map onto the set of consensus peaks
Columns specified here will be simplified into a single column. Should only be character or factor columns
String to separate values when collapsing columns
String to separate values when adding column names
logical(1) Include the original ranges as character columns
Passed to makeConsensus
GRanges object with a set of consensus ranges across all list elements and values from each element mapped to these consensus ranges.
If requested (include = TRUE
) the original ranges are returned as
character columns, as there will be multiple NA values in each.
Starting with a GRangesList, make a single GRanges object with select columns from each element added to the new object
a <- GRanges(paste0("chr1:", seq(1, 61, by = 20)))
width(a) <- 5
a$logFC <- rnorm(length(a))
a_g <- as.list(paste("Gene", seq_along(a)))
a_g[[1]] <- c("Gene 0", a_g[[1]])
a$genes <- as(a_g, "CompressedList")
b <- GRanges("chr1:61-70")
b$logFC <- rnorm(1)
b$genes <- as(list("Gene 5"), "CompressedList")
grl <- GRangesList(A = a, B = b)
mapGrlCols(grl, var = "logFC")
#> GRanges object with 4 ranges and 2 metadata columns:
#> seqnames ranges strand | A_logFC B_logFC
#> <Rle> <IRanges> <Rle> | <numeric> <numeric>
#> [1] chr1 1-5 * | -0.1777754 NA
#> [2] chr1 21-25 * | 0.0419013 NA
#> [3] chr1 41-45 * | -0.1453987 NA
#> [4] chr1 61-70 * | 0.6335019 0.839524
#> -------
#> seqinfo: 1 sequence from an unspecified genome; no seqlengths
## This forms a union of overlapping rangesby default
## Pass methods to makeConsensus() to change to regions with coverage == 2
mapGrlCols(grl, var = "logFC", method = "coverage", p = 1)
#> GRanges object with 1 range and 2 metadata columns:
#> seqnames ranges strand | A_logFC B_logFC
#> <Rle> <IRanges> <Rle> | <numeric> <numeric>
#> [1] chr1 61-65 * | 0.633502 0.839524
#> -------
#> seqinfo: 1 sequence from an unspecified genome
## Columns can be collapsed to merge into a single column
mapGrlCols(grl, var = "logFC", collapse = "genes")
#> GRanges object with 4 ranges and 3 metadata columns:
#> seqnames ranges strand | A_logFC B_logFC genes
#> <Rle> <IRanges> <Rle> | <numeric> <numeric> <character>
#> [1] chr1 1-5 * | -0.1777754 NA Gene 0; Gene 1
#> [2] chr1 21-25 * | 0.0419013 NA Gene 2
#> [3] chr1 41-45 * | -0.1453987 NA Gene 3
#> [4] chr1 61-70 * | 0.6335019 0.839524 Gene 4; Gene 5
#> -------
#> seqinfo: 1 sequence from an unspecified genome; no seqlengths
## Original ranges can also be included
mapGrlCols(grl, collapse = "genes", include = TRUE)
#> GRanges object with 4 ranges and 3 metadata columns:
#> seqnames ranges strand | A B genes
#> <Rle> <IRanges> <Rle> | <character> <character> <character>
#> [1] chr1 1-5 * | chr1:1-5 <NA> Gene 0; Gene 1
#> [2] chr1 21-25 * | chr1:21-25 <NA> Gene 2
#> [3] chr1 41-45 * | chr1:41-45 <NA> Gene 3
#> [4] chr1 61-70 * | chr1:61-65 chr1:61-70 Gene 4; Gene 5
#> -------
#> seqinfo: 1 sequence from an unspecified genome; no seqlengths