-
Notifications
You must be signed in to change notification settings - Fork 25
/
Copy pathflank.cpp
158 lines (123 loc) · 3.76 KB
/
flank.cpp
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
// flank.cpp
//
// Copyright (C) 2016 - 2018 Jay Hesselberth and Kent Riemondy
//
// This file is part of valr.
//
// This software may be modified and distributed under the terms
// of the MIT license. See the LICENSE file for details.
#include "valr.h"
void check_coords(int start, int end,
int chrom_size, int idx, bool trim,
std::vector<int>& starts_out,
std::vector<int>& ends_out,
std::vector<int>& df_idx) {
if (start == end) return ;
if (start >= 0 && end <= chrom_size) {
starts_out.push_back(start);
ends_out.push_back(end);
df_idx.push_back(idx);
} else if (trim) {
if (start < 0) {
starts_out.push_back(0) ;
} else {
starts_out.push_back(start) ;
}
if (end > chrom_size) {
ends_out.push_back(chrom_size) ;
} else {
ends_out.push_back(end) ;
}
df_idx.push_back(idx);
} // else trim
}
//[[Rcpp::export]]
DataFrame flank_impl(DataFrame df, DataFrame genome,
double both = 0, double left = 0, double right = 0,
bool fraction = false, bool stranded = false, bool trim = false) {
std::vector<std::string> chroms = df["chrom"];
IntegerVector starts = df["start"];
IntegerVector ends = df["end"];
// storage for outputs
std::vector<int> starts_out;
std::vector<int> ends_out;
std::vector<int> df_idx;
genome_map_t chrom_sizes = makeChromSizes(genome);
int lstart, lend, rstart, rend ;
if (stranded) {
std::vector<std::string> strands = df["strand"];
for (int i = 0; i < starts.size(); i++) {
int start = starts[i] ;
int end = ends[i] ;
double size = end - start;
if (fraction) {
if (strands[i] == "+") {
lstart = start - std::round(size * left);
lend = start;
rstart = end;
rend = end + std::round(size * right);
} else {
lstart = end;
lend = end + std::round(size * left) ;
rstart = start - std::round(size * right) ;
rend = start ;
}
} else {
if (strands[i] == "+") {
lstart = start - left;
lend = start;
rstart = end;
rend = end + right;
} else {
lstart = end;
lend = end + left;
rstart = start - right;
rend = start;
}
}
std::string chrom = chroms[i];
int chrom_size = chrom_sizes[chrom];
// check and save coordinates
check_coords(lstart, lend, chrom_size, i, trim,
starts_out, ends_out, df_idx) ;
check_coords(rstart, rend, chrom_size, i, trim,
starts_out, ends_out, df_idx) ;
}
} else { // no strand
for (int i = 0; i < starts.size(); i++) {
int start = starts[i] ;
int end = ends[i] ;
double size = end - start;
if (fraction) {
lstart = start - std::round(size * left);
lend = start;
rstart = end;
rend = end + std::round(size * right);
} else {
lstart = start - left;
lend = start;
rstart = end;
rend = end + right;
}
std::string chrom = chroms[i];
int chrom_size = chrom_sizes[chrom];
// check and save coordinates
check_coords(lstart, lend, chrom_size, i, trim,
starts_out, ends_out, df_idx) ;
check_coords(rstart, rend, chrom_size, i, trim,
starts_out, ends_out, df_idx) ;
}
}
DataFrame out = subset_dataframe(df, df_idx) ;
out["start"] = starts_out;
out["end"] = ends_out;
return out;
}
/*** R
library(valr)
library(dplyr)
genome <- read_genome(valr_example('hg19.chrom.sizes.gz'))
x <- bed_random(genome)
devtools::load_all()
flank_impl(x, genome, both = 100) %>% as_data_frame()
*/