forked from BioJulia/BioSequences.jl
-
Notifications
You must be signed in to change notification settings - Fork 0
Expand file tree
/
Copy pathBioSequences.jl
More file actions
312 lines (260 loc) · 5.79 KB
/
BioSequences.jl
File metadata and controls
312 lines (260 loc) · 5.79 KB
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
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
### BioSequences.jl
###
### A julia package for the representation and manipulation of biological sequences.
###
### This file is a part of BioJulia.
### License is MIT: https://github.com/BioJulia/BioSequences.jl/blob/master/LICENSE
module BioSequences
export
###
### Symbols
###
# Types & aliases
NucleicAcid,
DNA,
RNA,
DNA_A,
DNA_C,
DNA_G,
DNA_T,
DNA_M,
DNA_R,
DNA_W,
DNA_S,
DNA_Y,
DNA_K,
DNA_V,
DNA_H,
DNA_D,
DNA_B,
DNA_N,
DNA_Gap,
ACGT,
ACGTN,
RNA_A,
RNA_C,
RNA_G,
RNA_U,
RNA_M,
RNA_R,
RNA_W,
RNA_S,
RNA_Y,
RNA_K,
RNA_V,
RNA_H,
RNA_D,
RNA_B,
RNA_N,
RNA_Gap,
ACGU,
ACGUN,
AminoAcid,
AA_A,
AA_R,
AA_N,
AA_D,
AA_C,
AA_Q,
AA_E,
AA_G,
AA_H,
AA_I,
AA_L,
AA_K,
AA_M,
AA_F,
AA_P,
AA_S,
AA_T,
AA_W,
AA_Y,
AA_V,
AA_O,
AA_U,
AA_B,
AA_J,
AA_Z,
AA_X,
AA_Term,
AA_Gap,
# Predicates
isGC,
iscompatible,
isambiguous,
iscertain,
isgap,
ispurine,
ispyrimidine,
###
### Alphabets
###
# Types & aliases
Alphabet,
NucleicAcidAlphabet,
DNAAlphabet,
RNAAlphabet,
AminoAcidAlphabet,
###
### BioSequences
###
join!,
# Type & aliases
BioSequence,
NucSeq,
AASeq,
# Predicates
ispalindromic,
hasambiguity,
isrepetitive,
iscanonical,
# Transformations
canonical,
canonical!,
complement,
complement!,
reverse_complement,
reverse_complement!,
ungap,
ungap!,
join!,
###
### LongSequence
###
# Type & aliases
LongSequence,
LongNuc,
LongDNA,
LongRNA,
LongAA,
LongSubSeq,
# Random
SamplerUniform,
SamplerWeighted,
randseq,
randdnaseq,
randrnaseq,
randaaseq,
###
### Sequence literals
###
@dna_str,
@rna_str,
@aa_str,
@biore_str,
@prosite_str,
BioRegex,
BioRegexMatch,
matched,
captured,
alphabet, # TODO: Resolve the use of alphabet - it's from BioSymbols.jl
symbols,
gap,
mismatches,
matches,
n_ambiguous,
n_gaps,
n_certain,
gc_content,
translate!,
translate,
ncbi_trans_table,
# Search
ExactSearchQuery,
ApproximateSearchQuery,
PFM,
PWM,
PWMSearchQuery,
maxscore,
scoreat,
seqmatrix,
majorityvote
using BioSymbols
import Twiddle: enumerate_nibbles,
count_0000_nibbles,
count_1111_nibbles,
count_nonzero_nibbles,
count_00_bitpairs,
count_01_bitpairs,
count_10_bitpairs,
count_11_bitpairs,
count_nonzero_bitpairs,
repeatpattern
using Random
BioSymbols.gap(::Type{Char}) = '-'
include("alphabet.jl")
# Load the bit-twiddling internals that optimised BioSequences methods depend on.
include("bit-manipulation/bit-manipulation.jl")
# The generic, abstract BioSequence type
include("biosequence/biosequence.jl")
# The definition of the LongSequence concrete type, and its method overloads...
include("longsequences/longsequence.jl")
include("longsequences/hash.jl")
include("longsequences/randseq.jl")
include("geneticcode.jl")
# Pattern searching in sequences...
include("search/ExactSearchQuery.jl")
include("search/ApproxSearchQuery.jl")
include("search/re.jl")
include("search/pwm.jl")
struct Search{Q,I}
query::Q
itr::I
overlap::Bool
end
const DEFAULT_OVERLAP = true
search(query, itr; overlap = DEFAULT_OVERLAP) = Search(query, itr, overlap)
function Base.iterate(itr::Search, state=firstindex(itr.itr))
val = findnext(itr.query, itr.itr, state)
val === nothing && return nothing
state = itr.overlap ? first(val) + 1 : last(val) + 1
return val, state
end
const HasRangeEltype = Union{<:ExactSearchQuery, <:ApproximateSearchQuery, <:Regex}
Base.eltype(::Type{<:Search{Q}}) where {Q<:HasRangeEltype} = UnitRange{Int}
Base.eltype(::Type{<:Search}) = Int
Base.IteratorSize(::Type{<:Search}) = Base.SizeUnknown()
"""
findall(pattern, sequence::BioSequence[,rng::UnitRange{Int}]; overlap::Bool=true)::Vector
Find all occurrences of `pattern` in `sequence`.
The return value is a vector of ranges of indices where the matching sequences were found.
If there are no matching sequences, the return value is an empty vector.
The search is restricted to the specified range when `rng` is set.
With the keyword argument `overlap` set as `true`, the start index for the next search gets set to the start of the current match plus one; if set to `false`, the start index for the next search gets set to the end of the current match plus one.
The default value for the keyword argument `overlap` is `true`.
The `pattern` can be a `Biosymbol` or a search query.
See also [`ExactSearchQuery`](@ref), [`ApproximateSearchQuery`](@ref), [`PWMSearchQuery`](@ref).
# Examples
```jldoctest
julia> seq = dna"ACACACAC";
julia> findall(DNA_A, seq)
4-element Vector{Int64}:
1
3
5
7
julia> findall(ExactSearchQuery(dna"ACAC"), seq)
3-element Vector{UnitRange{Int64}}:
1:4
3:6
5:8
julia> findall(ExactSearchQuery(dna"ACAC"), seq; overlap=false)
2-element Vector{UnitRange{Int64}}:
1:4
5:8
julia> findall(ExactSearchQuery(dna"ACAC"), seq, 2:7; overlap=false)
1-element Vector{UnitRange{Int64}}:
3:6
```
"""
function Base.findall(pat, seq::BioSequence; overlap::Bool = DEFAULT_OVERLAP)
return collect(search(pat, seq; overlap))
end
# Fix ambiguity with Base's findall
Base.findall(f::Function, seq::BioSequence) = collect(search(f, seq; overlap=DEFAULT_OVERLAP))
function Base.findall(pat, seq::BioSequence, rng::UnitRange{Int}; overlap::Bool = DEFAULT_OVERLAP)
v = view(seq, rng)
itr = search(pat, v; overlap)
return map(x->parentindices(v)[1][x], itr)
end
end # module BioSequences