forked from jimhester/GenomicRanges
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathcoverage-methods.R
86 lines (78 loc) · 3.12 KB
/
coverage-methods.R
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
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
### =========================================================================
### "coverage" methods
### -------------------------------------------------------------------------
###
### TODO: Merge with Biostrings:::.V_recycle() and put in IRanges.
.recycle <- function(x, skeleton_len, x.label, skeleton.label)
{
x_len <- length(x)
if (x_len == skeleton_len)
return(x)
if (x_len < skeleton_len) {
if (x_len == 0L)
stop("cannot recycle zero-length '", x.label, "' ",
"to the length of '", skeleton.label, "'")
} else {
if (x_len >= 2L)
stop("'", x.label, "' is longer than '", skeleton.label, "'")
}
if (skeleton_len %% x_len != 0L)
warning("'", x.label, "' length is not a divisor ",
"of '", skeleton.label, "' length")
rep(x, length.out=skeleton_len)
}
### Returns a list-like object.
.normarg_shift_or_weight <- function(arg, arg.label, x)
{
if (is.list(arg) || is(arg, "List")) {
if (!identical(names(arg), seqlevels(x)))
stop("when '", arg.label, "' is a list-like object, it must ",
"have 1 list element per seqlevel in 'x', and its names ",
"must be exactly 'seqlevels(x)'")
return(arg)
}
if (isSingleString(arg)) {
x_mcols <- mcols(x)
if (!is(x_mcols, "DataTable")
|| sum(colnames(x_mcols) == arg) != 1L)
stop("'mcols(x)' has 0 or more than 1 \"",
arg, "\" columns")
arg <- x_mcols[ , arg]
}
if (!(is.numeric(arg) || is(arg, "Rle") && is.numeric(runValue(arg))))
stop("'", arg.label, "' must be a numeric vector, a single string, ",
"or a list-like object")
split(.recycle(arg, length(x), arg.label, "x"), seqnames(x))
}
setMethod("coverage", "GenomicRanges",
function(x, shift=0L, width=NULL, weight=1L,
method=c("auto", "sort", "hash"))
{
## Normalize 'shift'.
shift <- .normarg_shift_or_weight(shift, "shift", x)
## Just handle the default 'width' here. Non default will be checked
## in IRanges:::.CompressedIRangesList.coverage().
if (is.null(width))
width <- seqlengths(x)
## Normalize 'weight'.
weight <- .normarg_shift_or_weight(weight, "weight", x)
x_ranges_list <- split(ranges(x), seqnames(x))
circle.length <- seqlengths(x)
circle.length[!(isCircular(x) %in% TRUE)] <- NA_integer_
IRanges:::.CompressedIRangesList.coverage(x_ranges_list,
shift=shift,
width=width,
weight=weight,
circle.length=circle.length,
method=method,
x_names.label="'seqlevels(x)'")
}
)
setMethod("coverage", "GRangesList",
function(x, shift=0L, width=NULL, weight=1L,
method=c("auto", "sort", "hash"))
{
coverage(x@unlistData, shift=shift, width=width, weight=weight,
method=method)
}
)