-
Notifications
You must be signed in to change notification settings - Fork 19
/
ranges-join-nearest.R
225 lines (200 loc) · 8.02 KB
/
ranges-join-nearest.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
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
#' Find nearest neighbours between two Ranges objects
#'
#' @param x,y Ranges objects, add the nearest neighbours of ranges in x to
#' those in y.
#' @param suffix A character vector of length two used to identify metadata columns
#' @param distance logical vector whether to add a column named "distance"
#' containing the distance to the nearest region. If set to a character vector
#' of length 1, will use that as distance column name.
#'
#' @details By default `join_nearest` will find arbitrary nearest
#' neighbours in either direction and ignore any strand information.
#' The `join_nearest_left` and `join_nearest_right` methods
#' will find arbitrary nearest neighbour ranges on x that are left/right of
#' those on y and ignore any strand information.
#'
#' The `join_nearest_upstream` method will find arbitrary nearest
#' neighbour ranges on x that are upstream of those on y. This takes into
#' account strandedness of the ranges.
#' On the positive strand nearest upstream will be on the
#' left and on the negative strand nearest upstream will be on the right.
#'
#' The `join_nearest_downstream` method will find arbitrary nearest
#' neighbour ranges on x that are upstream of those on y. This takes into
#' account strandedness of the ranges.On the positive strand nearest downstream
#' will be on the right and on the negative strand nearest upstream will be on
#' the left.
#'
#' @return A Ranges object corresponding to the nearest ranges, all metadata
#' is copied over from the right-hand side ranges `y`.
#'
#' @examples
#' query <- data.frame(start = c(5,10, 15,20),
#' width = 5,
#' gc = runif(4)) %>%
#' as_iranges()
#' subject <- data.frame(start = c(2:6, 24),
#' width = 3:8,
#' label = letters[1:6]) %>%
#' as_iranges()
#'
#' join_nearest(query, subject)
#' join_nearest_left(query, subject)
#' join_nearest_right(query, subject)
#'
#' subject <- data.frame(seqnames = "chr1",
#' start = c(11,101),
#' end = c(21, 200),
#' name = c("a1", "a2"),
#' strand = c("+", "-"),
#' score = c(1,2)) %>%
#' as_granges()
#' query <- data.frame(seqnames = "chr1",
#' strand = c("+", "-", "+", "-"),
#' start = c(21,91,101,201),
#' end = c(30,101,110,210),
#' name = paste0("b", 1:4),
#' score = 1:4) %>%
#' as_granges()
#' join_nearest_upstream(query, subject)
#' join_nearest_downstream(query, subject)
#' @rdname ranges-nearest
#' @importFrom IRanges nearest
#' @export
join_nearest <- function(x,y, suffix = c(".x", ".y"), distance = FALSE) {
hits <- hits_nearest(x, y)
distance_col <- set_distance_col(distance)
expand_by_hits(x, y, suffix, hits, hits_mcols_to_keep = distance_col)
}
#' @rdname ranges-nearest
#' @export
join_nearest_left <- function(x,y, suffix = c(".x", ".y"), distance = FALSE) {
hits <- hits_nearest_left(x, y)
distance_col <- set_distance_col(distance)
expand_by_hits(x, y, suffix, hits, hits_mcols_to_keep = distance_col)
}
#' @rdname ranges-nearest
#' @export
join_nearest_right <- function(x, y, suffix = c(".x", ".y"), distance = FALSE) {
hits <- hits_nearest_right(x, y)
distance_col <- set_distance_col(distance)
expand_by_hits(x, y, suffix, hits, hits_mcols_to_keep = distance_col)
}
#' @rdname ranges-nearest
#' @export
join_nearest_upstream <- function(x, y, suffix = c(".x", ".y"), distance = FALSE) { UseMethod("join_nearest_upstream")}
#' @export
join_nearest_upstream.GenomicRanges <- function(x, y, suffix = c(".x", ".y"), distance = FALSE) {
hits <- hits_nearest_upstream(x, y)
distance_col <- set_distance_col(distance)
expand_by_hits(x, y, suffix, hits, hits_mcols_to_keep = distance_col)
}
#' @rdname ranges-nearest
#' @importFrom IRanges nearest
#' @export
join_nearest_downstream <- function(x, y, suffix = c(".x", ".y"), distance = FALSE) { UseMethod("join_nearest_downstream")}
#' @export
join_nearest_downstream.GenomicRanges <- function(x, y, suffix = c(".x", ".y"), distance = FALSE) {
hits <- hits_nearest_downstream(x, y)
distance_col <- set_distance_col(distance)
expand_by_hits(x, y, suffix, hits, hits_mcols_to_keep = distance_col)
}
# -----------------------
# Directed hits helpers
hits_nearest <- function(x, y){
UseMethod("hits_nearest")
}
hits_nearest_left <- function(x, y){
UseMethod("hits_nearest_left")
}
hits_nearest_right <- function(x, y){
UseMethod("hits_nearest_right")
}
hits_nearest_upstream <- function(x, y){
UseMethod("hits_nearest_upstream")
}
hits_nearest_downstream <- function(x, y){
UseMethod("hits_nearest_downstream")
}
#' @importFrom IRanges distanceToNearest
hits_nearest.IRanges <- function(x, y){
make_hits(x, y, distanceToNearest, select = "arbitrary")
}
#' @importFrom GenomicRanges distanceToNearest
hits_nearest.GenomicRanges <- function(x, y){
make_hits(x, y, distanceToNearest, select = "arbitrary", ignore.strand = TRUE)
}
#' @importFrom IRanges distanceToNearest
hits_nearest_left.IRanges <- function(x, y){
hits <- make_hits(x, y, distanceToNearest, select = "all")
get_hits_left(x, y, hits)
}
#' @importFrom GenomicRanges distanceToNearest
hits_nearest_left.GenomicRanges <- function(x, y){
hits <- make_hits(x, y, distanceToNearest, select = "all", ignore.strand = TRUE)
get_hits_left(x,y, hits)
}
#' @importFrom IRanges distanceToNearest
hits_nearest_right.IRanges <- function(x, y){
hits <- make_hits(x, y, distanceToNearest, select = "all")
get_hits_right(x, y, hits)
}
#' @importFrom GenomicRanges distanceToNearest
hits_nearest_right.GenomicRanges <- function(x, y){
hits <- make_hits(x, y, distanceToNearest, select = "all", ignore.strand = TRUE)
get_hits_right(x, y, hits)
}
#' @importFrom GenomicRanges distanceToNearest
hits_nearest_upstream.GenomicRanges <- function(x, y){
hits <- distanceToNearest(x,y, select = "all", ignore.strand = FALSE)
get_hits_upstream(x, y, hits)
}
#' @importFrom GenomicRanges distanceToNearest
hits_nearest_downstream.GenomicRanges <- function(x, y){
hits <- distanceToNearest(x,y, select = "all", ignore.strand = FALSE)
get_hits_downstream(x, y, hits)
}
# -----------------------
# Hits logic helpers
get_hits_left <- function(x, y, hits){
mcols(hits)$is_left <- start(y[subjectHits(hits)]) <= start(x[queryHits(hits)]) &
end(y[subjectHits(hits)]) <= start(x[queryHits(hits)])
hits[mcols(hits)$is_left]
}
get_hits_right <- function(x, y, hits){
mcols(hits)$is_right <- end(x[queryHits(hits)]) <= start(y[subjectHits(hits)])
hits[mcols(hits)$is_right]
}
get_hits_upstream <- function(x, y, hits){
mcols(hits)$is_right <- end(x[queryHits(hits)]) <= start(y[subjectHits(hits)])
mcols(hits)$is_left <- start(x[queryHits(hits)]) >= start(y[subjectHits(hits)])
mcols(hits)$direction <- strand(x[queryHits(hits)])
mcols(hits)$is_upstream <- ifelse(mcols(hits)$direction == "+",
mcols(hits)$is_left,
mcols(hits)$is_right)
hits[mcols(hits)$is_upstream]
}
get_hits_downstream <- function(x, y, hits){
mcols(hits)$is_right <- end(x[queryHits(hits)]) <= start(y[subjectHits(hits)])
mcols(hits)$is_left <- start(x[queryHits(hits)]) >= start(y[subjectHits(hits)])
mcols(hits)$direction <- strand(x[queryHits(hits)])
# on positive strand nearest downtream will be on the right
# on negative strand nearest downstream will be on the left
mcols(hits)$is_downstream <- ifelse(mcols(hits)$direction == "-",
mcols(hits)$is_left,
mcols(hits)$is_right)
hits[mcols(hits)$is_downstream]
}
set_distance_col <- function(distance){
if (distance == TRUE){
keep_distance <- "distance"
} else if (is.character(distance) & length(distance) == 1) {
keep_distance <- c("distance")
names(keep_distance) <- distance
} else if (distance == FALSE) {
keep_distance = NULL
} else {
stop("Invalid argument passed to distance")
}
return(keep_distance)
}