mergeGRangesData: Merge GRanges objects

Description Usage Arguments Value Subsetting a multiplexed GRanges object Author(s) See Also Examples

View source: R/dataset_functions.R

Description

Merges 2 or more GRanges objects by combining all of their ranges and associated signal (e.g. readcounts). If multiplex = TRUE, the input datasets are reversibly combined into a multiplexed GRanges containing a field for each input dataset.

Usage

1
2
3
4
5
6
7
8
mergeGRangesData(
  ...,
  field = "score",
  multiplex = FALSE,
  makeBRG = TRUE,
  exact_overlaps = FALSE,
  ncores = getOption("mc.cores", 2L)
)

Arguments

...

Any number of GRanges objects in which signal (e.g. readcounts) are contained within metadata. Lists of GRanges can also be passed, but they must be named lists if multiplex = TRUE. Multiple lists can be passed, but if any inputs are lists, then all inputs must be lists.

field

One or more input metadata fields to be combined, typically the "score" field. Fields typically contain coverage information. If only a single field is given (i.e. all input GRanges use the same field), that same field will be used for the output. Otherwise, the "score" metadata field will be used by default. The output metadata fields are different if multiplex is enabled.

multiplex

When set to FALSE (the default), input GRanges are merged irreversibly into a single new GRange, effectively combining the reads from different experiments. When multiplex = TRUE, the input GRanges data are reversibly combined into a multiplexed GRanges object, such that each input GRanges object has its own metadata field in the output.

makeBRG

If TRUE (the default), the output GRanges will be made "basepair-resolution" via makeGRangesBRG. This is always set to codeFALSE if exact_overlaps = TRUE.

exact_overlaps

By default (FALSE), any ranges that cover more than a single base are treated as "coverage" signal (see getCountsByRegions). If exact_overlaps = TRUE, all input ranges are conserved exactly as they are; ranges are only combined if they overlap exactly, and the signal of any combined range is the sum of the input ranges that were merged.

ncores

Number of cores to use for computations.

Value

A disjoint, basepair-resolution (single-width) GRanges object comprised of all ranges found in the input GRanges objects.

If multiplex = FALSE, single fields from each input are combined into a single field in the output, the total signal of which is the sum of all input GRanges.

If multiplex = TRUE, each field of the output corresponds to an input GRanges object.

Subsetting a multiplexed GRanges object

If multiplex = TRUE, the datasets are only combined into a single object, but the data themselves are not combined. To subset field_i, corresponding to input dataset_i:

multi.gr <- mergeGRangesData(gr1, gr2, multiplex = TRUE) subset(multi.gr, gr1 != 0, select = gr1) # select gr1

Author(s)

Mike DeBerardine

See Also

makeGRangesBRG

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
data("PROseq") # load included PROseq data

#--------------------------------------------------#
# divide & recombine PROseq (no overlapping positions)
#--------------------------------------------------#

thirds <- floor( (1:3)/3 * length(PROseq) )
ps_1 <- PROseq[1:thirds[1]]
ps_2 <- PROseq[(thirds[1]+1):thirds[2]]
ps_3 <- PROseq[(thirds[2]+1):thirds[3]]

# re-merge
length(PROseq)
length(ps_1)
length(mergeGRangesData(ps_1, ps_2, ncores = 1))
length(mergeGRangesData(ps_1, ps_2, ps_3, ncores = 1))

#--------------------------------------------------#
# combine PRO-seq with overlapping positions
#--------------------------------------------------#

gr1 <- PROseq[10:13]
gr2 <- PROseq[12:15]

PROseq[10:15]

mergeGRangesData(gr1, gr2, ncores = 1)

#--------------------------------------------------#
# multiplex separate PRO-seq experiments
#--------------------------------------------------#

multi.gr <- mergeGRangesData(gr1, gr2, multiplex = TRUE, ncores = 1)
multi.gr

#--------------------------------------------------#
# subset a multiplexed GRanges object
#--------------------------------------------------#

subset(multi.gr, gr1 > 0)

subset(multi.gr, gr1 > 0, select = gr1)

BRGenomics documentation built on Nov. 8, 2020, 8:03 p.m.