0
|
1 #!/usr/bin/env python3
|
|
2
|
|
3 import numpy as np
|
|
4 import subprocess
|
|
5 import math
|
|
6 import time
|
|
7 from operator import itemgetter
|
|
8 from collections import Counter
|
|
9 from itertools import groupby
|
|
10 import os
|
|
11 import configuration
|
|
12 from tempfile import NamedTemporaryFile
|
|
13 import sys
|
|
14 import warnings
|
|
15 import shutil
|
|
16 from collections import defaultdict
|
|
17
|
|
18 np.set_printoptions(threshold=sys.maxsize)
|
|
19
|
|
20
|
|
21 def alignment_scoring():
|
|
22 ''' Create hash table for alignment similarity counting: for every
|
|
23 combination of aminoacids in alignment assign score from protein
|
|
24 scoring matrix defined in configuration file '''
|
|
25 score_dict = {}
|
|
26 with open(configuration.SC_MATRIX) as smatrix:
|
|
27 count = 1
|
|
28 for line in smatrix:
|
|
29 if not line.startswith("#"):
|
|
30 if count == 1:
|
|
31 aa_all = line.rstrip().replace(" ", "")
|
|
32 else:
|
|
33 count_aa = 1
|
|
34 line = list(filter(None, line.rstrip().split(" ")))
|
|
35 for aa in aa_all:
|
|
36 score_dict["{}{}".format(line[0], aa)] = line[count_aa]
|
|
37 count_aa += 1
|
|
38 count += 1
|
|
39 return score_dict
|
|
40
|
|
41
|
|
42 def characterize_fasta(QUERY, WIN_DOM):
|
|
43 ''' Find the sequences, their lengths, starts, ends and if
|
|
44 they exceed the window '''
|
|
45 with open(QUERY) as query:
|
|
46 headers = []
|
|
47 fasta_lengths = []
|
|
48 seq_starts = []
|
|
49 seq_ends = []
|
|
50 fasta_chunk_len = 0
|
|
51 count_line = 1
|
|
52 for line in query:
|
|
53 line = line.rstrip()
|
|
54 if line.startswith(">"):
|
|
55 headers.append(line.rstrip())
|
|
56 fasta_lengths.append(fasta_chunk_len)
|
|
57 fasta_chunk_len = 0
|
|
58 seq_starts.append(count_line + 1)
|
|
59 seq_ends.append(count_line - 1)
|
|
60 else:
|
|
61 fasta_chunk_len += len(line)
|
|
62 count_line += 1
|
|
63 seq_ends.append(count_line)
|
|
64 seq_ends = seq_ends[1:]
|
|
65 fasta_lengths.append(fasta_chunk_len)
|
|
66 fasta_lengths = fasta_lengths[1:]
|
|
67 # control if there are correct (unique) names for individual seqs:
|
|
68 # LASTAL takes seqs IDs till the first space which can then create problems with ambiguous records
|
|
69 if len(headers) > len(set([header.split(" ")[0] for header in headers
|
|
70 ])):
|
|
71 raise NameError(
|
|
72 '''Sequences in multifasta format are not named correctly:
|
|
73 seq IDs (before the first space) are the same''')
|
|
74
|
|
75 above_win = [idx
|
|
76 for idx, value in enumerate(fasta_lengths) if value > WIN_DOM]
|
|
77 below_win = [idx
|
|
78 for idx, value in enumerate(fasta_lengths)
|
|
79 if value <= WIN_DOM]
|
|
80 lens_above_win = np.array(fasta_lengths)[above_win]
|
|
81 return headers, above_win, below_win, lens_above_win, seq_starts, seq_ends
|
|
82
|
|
83
|
|
84 def split_fasta(QUERY, WIN_DOM, step, headers, above_win, below_win,
|
|
85 lens_above_win, seq_starts, seq_ends):
|
|
86 ''' Create temporary file containing all sequences - the ones that exceed
|
|
87 the window are cut with a set overlap (greater than domain size with a reserve) '''
|
|
88 with open(QUERY, "r") as query:
|
|
89 count_fasta_divided = 0
|
|
90 count_fasta_not_divided = 0
|
|
91 ntf = NamedTemporaryFile(delete=False)
|
|
92 divided = np.array(headers)[above_win]
|
|
93 row_length = configuration.FASTA_LINE
|
|
94 for line in query:
|
|
95 line = line.rstrip()
|
|
96 if line.startswith(">") and line in divided:
|
|
97 stop_line = seq_ends[above_win[
|
|
98 count_fasta_divided]] - seq_starts[above_win[
|
|
99 count_fasta_divided]] + 1
|
|
100 count_line = 0
|
|
101 whole_seq = []
|
|
102 for line2 in query:
|
|
103 whole_seq.append(line2.rstrip())
|
|
104 count_line += 1
|
|
105 if count_line == stop_line:
|
|
106 break
|
|
107 whole_seq = "".join(whole_seq)
|
|
108 ## create list of starting positions for individual parts of a seq with a step given by a window and overlap
|
|
109 windows_starts = list(range(0, lens_above_win[
|
|
110 count_fasta_divided], step))
|
|
111 ## create list of ending positions (starting pos + window), the last element is the whole seq length
|
|
112 windows_ends = [
|
|
113 x + WIN_DOM
|
|
114 if x + WIN_DOM < lens_above_win[count_fasta_divided] else
|
|
115 lens_above_win[count_fasta_divided] for x in windows_starts
|
|
116 ]
|
|
117 count_part = 1
|
|
118 for start_part, end_part in zip(windows_starts, windows_ends):
|
|
119 seq_part = whole_seq[start_part:end_part]
|
|
120 if count_part == len(windows_starts):
|
|
121 ntf.write("{}_DANTE_PART{}_LAST:{}-{}\n{}\n".format(
|
|
122 line.split(" ")[0], count_part, start_part + 1,
|
|
123 end_part, "\n".join([seq_part[i:i + row_length]
|
|
124 for i in range(0, len(
|
|
125 seq_part), row_length)
|
|
126 ])).encode("utf-8"))
|
|
127 else:
|
|
128 ntf.write("{}_DANTE_PART{}:{}-{}\n{}\n".format(
|
|
129 line.split(" ")[0], count_part, start_part + 1,
|
|
130 end_part, "\n".join([seq_part[i:i + row_length]
|
|
131 for i in range(0, len(
|
|
132 seq_part), row_length)
|
|
133 ])).encode("utf-8"))
|
|
134 count_part += 1
|
|
135 count_fasta_divided += 1
|
|
136 elif line.startswith(">") and line not in divided:
|
|
137 length_seq = seq_ends[below_win[
|
|
138 count_fasta_not_divided]] - seq_starts[below_win[
|
|
139 count_fasta_not_divided]] + 1
|
|
140 ntf.write("{}\n{}".format(line, "".join([query.readline(
|
|
141 ) for x in range(length_seq)])).encode("utf-8"))
|
|
142 count_fasta_not_divided += 1
|
|
143 query_temp = ntf.name
|
|
144 ntf.close()
|
|
145 return query_temp
|
|
146
|
|
147
|
|
148 def domain_annotation(elements, CLASSIFICATION):
|
|
149 ''' Assign protein domain to each hit from protein database '''
|
|
150 domains = []
|
|
151 annotations = []
|
|
152 with open(CLASSIFICATION, "r") as cl_tbl:
|
|
153 annotation = {}
|
|
154 for line in cl_tbl:
|
|
155 record = line.rstrip().split("\t")
|
|
156 annotation[record[0]] = record[1:]
|
|
157 for i in range(len(elements)):
|
|
158 domains.append(elements[i].split("__")[0].split("-")[1])
|
|
159 element_name = "__".join(elements[i].split("__")[1:])
|
|
160 if element_name in annotation.keys():
|
|
161 annotations.append("|".join([elements[i].split("__")[0].split("-")[
|
|
162 1], ("|".join(annotation[element_name]))]))
|
|
163 else:
|
|
164 annotations.append("unknown|unknown")
|
|
165 return annotations
|
|
166
|
|
167
|
|
168 def hits_processing(seq_len, start, end, strand):
|
|
169 ''' Gain hits intervals separately for forward and reverse strand '''
|
|
170 reverse_strand_idx = np.where(strand == "-")[0]
|
|
171 if not reverse_strand_idx.any():
|
|
172 start_pos_plus = start + 1
|
|
173 end_pos_plus = end
|
|
174 regions_plus = list(zip(start_pos_plus, end_pos_plus))
|
|
175 regions_minus = []
|
|
176 else:
|
|
177 reverse_strand_idx = reverse_strand_idx[0]
|
|
178 start_pos_plus = start[0:reverse_strand_idx] + 1
|
|
179 end_pos_plus = end[0:reverse_strand_idx]
|
|
180 start_pos_minus = seq_len[0] - end[reverse_strand_idx:] + 1
|
|
181 end_pos_minus = seq_len[0] - start[reverse_strand_idx:]
|
|
182 regions_plus = list(zip(start_pos_plus, end_pos_plus))
|
|
183 regions_minus = list(zip(start_pos_minus, end_pos_minus))
|
|
184 return reverse_strand_idx, regions_plus, regions_minus
|
|
185
|
|
186
|
|
187 def overlapping_regions(input_data):
|
|
188 ''' Join all overlapping intervals(hits) to clusters (potential domains),
|
|
189 get list of start-end positions of individual hits within the interval,
|
|
190 list of minimus and maximums as well as the indices in the original
|
|
191 sequence_hits structure for the hits belonging to the same clusters '''
|
|
192 if input_data:
|
|
193 sorted_idx, sorted_data = zip(*sorted(
|
|
194 [(index, data) for index, data in enumerate(input_data)],
|
|
195 key=itemgetter(1)))
|
|
196 merged_ends = input_data[sorted_idx[0]][1]
|
|
197 intervals = []
|
|
198 data = []
|
|
199 output_intervals = []
|
|
200 output_data = []
|
|
201 for i, j in zip(sorted_idx, sorted_data):
|
|
202 if input_data[i][0] < merged_ends:
|
|
203 merged_ends = max(input_data[i][1], merged_ends)
|
|
204 intervals.append(i)
|
|
205 data.append(j)
|
|
206 else:
|
|
207 output_intervals.append(intervals)
|
|
208 output_data.append(data)
|
|
209 intervals = []
|
|
210 data = []
|
|
211 intervals.append(i)
|
|
212 data.append(j)
|
|
213 merged_ends = input_data[i][1]
|
|
214 output_intervals.append(intervals)
|
|
215 output_data.append(data)
|
|
216 mins = [x[0][0] for x in output_data]
|
|
217 maxs = [max(x, key=itemgetter(1))[1] for x in output_data]
|
|
218 else:
|
|
219 mins = []
|
|
220 maxs = []
|
|
221 output_intervals = []
|
|
222 output_data = []
|
|
223 return mins, maxs, output_data, output_intervals
|
|
224
|
|
225
|
|
226 def annotations_dict(annotations):
|
|
227 ''' Hash table where annotations of the hits within a clusters are the keys.
|
|
228 Each annotation has serial number assigned which indexes the row in the score_table '''
|
|
229 classes_dict = {classes: idx
|
|
230 for idx, classes in enumerate(set(annotations))}
|
|
231 return classes_dict
|
|
232
|
|
233
|
|
234 def score_table(mins, maxs, data, annotations, scores, CLASSIFICATION):
|
|
235 ''' Score table is created based on the annotations occurance in the cluster.
|
|
236 Matrix axis y corresponds to individual annotations (indexed according to classes_dict),
|
|
237 axis x represents positions of analyzed seq in a given cluster.
|
|
238 For every hit within cluster, array of scores on the corresponding position
|
|
239 is recorded to the table in case if the score on certain position is so far the highest
|
|
240 for the certain position and certain annotation '''
|
|
241 classes_dict = annotations_dict(annotations)
|
|
242 score_matrix = np.zeros((len(classes_dict), maxs - mins + 1), dtype=int)
|
|
243 count = 0
|
|
244 for item in annotations:
|
|
245 saved_scores = score_matrix[classes_dict[item], data[count][0] - mins:
|
|
246 data[count][1] - mins + 1]
|
|
247 new_scores = [scores[count]] * len(saved_scores)
|
|
248 score_matrix[classes_dict[item], data[count][0] - mins:data[count][
|
|
249 1] - mins + 1] = [max(*pos_score)
|
|
250 for pos_score in zip(saved_scores, new_scores)]
|
|
251 count += 1
|
|
252 return score_matrix, classes_dict
|
|
253
|
|
254
|
|
255 def score_matrix_evaluation(score_matrix, classes_dict, THRESHOLD_SCORE):
|
|
256 ''' Score matrix is evaluated based on each position.
|
|
257 For every position the list of annotations with a score which reaches
|
|
258 certain percentage of the overal best score of the cluster are stored '''
|
|
259 ann_per_reg = []
|
|
260 overal_best_score_reg = max((score_matrix.max(axis=1)))
|
|
261 for position in score_matrix.T:
|
|
262 ## score threshold calculated as a percentage of the OVERALL best score in the cluster
|
|
263 threshold = overal_best_score_reg * THRESHOLD_SCORE / 100
|
|
264 above_th = [idx
|
|
265 for idx, score in enumerate(position)
|
|
266 if position[idx] >= threshold]
|
|
267 ## select unique annotations in one position that are above threshold
|
|
268 ann_per_pos = list(set(
|
|
269 [key for key, value in classes_dict.items() if value in above_th]))
|
|
270 ann_per_reg.append(ann_per_pos)
|
|
271 return ann_per_reg
|
|
272
|
|
273
|
|
274 def group_annot_regs(ann_per_reg):
|
|
275 ''' Get list of domains, annotations, longest common annotations and
|
|
276 counts of positions with certain annotation per regions '''
|
|
277 ## tranform list of lists (potential multiple annotations for every position ) to flat list of all annotations
|
|
278 all_annotations = [item for sublist in ann_per_reg for item in sublist]
|
|
279 unique_annotations = list(set(all_annotations))
|
|
280 ann_pos_counts = [all_annotations.count(x) for x in unique_annotations]
|
|
281 unique_annotations = list(set(
|
|
282 [item for sublist in ann_per_reg for item in sublist]))
|
|
283 domain_type = list(set([annotation.split("|")[0]
|
|
284 for annotation in unique_annotations]))
|
|
285 classification_list = [annotation.split("|")
|
|
286 for annotation in unique_annotations]
|
|
287 ann_substring = "|".join(os.path.commonprefix(classification_list))
|
|
288 domain_type = "/".join(domain_type)
|
|
289 return domain_type, ann_substring, unique_annotations, ann_pos_counts
|
|
290
|
|
291
|
|
292 def best_score(scores, region):
|
|
293 ''' From overlapping intervals take the one with the highest score '''
|
|
294 ## if more hits have the same best score take only the first one
|
|
295 best_idx = region[np.where(scores == max(scores))[0][0]]
|
|
296 best_idx_reg = np.where(scores == max(scores))[0][0]
|
|
297 return best_idx, best_idx_reg
|
|
298
|
|
299
|
|
300 def create_gff3(domain_type, ann_substring, unique_annotations, ann_pos_counts,
|
|
301 dom_start, dom_end, step, best_idx, annotation_best,
|
|
302 db_name_best, db_starts_best, db_ends_best, strand, score,
|
|
303 seq_id, db_seq, query_seq, domain_size, positions, gff):
|
|
304 ''' Record obtained information about domain corresponding to individual cluster to common gff file '''
|
|
305 best_start = positions[best_idx][0]
|
|
306 best_end = positions[best_idx][1]
|
|
307 best_score = score[best_idx]
|
|
308 ## proportion of length of the best hit to the whole region length found by base
|
|
309 length_proportion = int((best_end - best_start + 1) /
|
|
310 (dom_end - dom_start + 1) * 100)
|
|
311 db_seq_best = db_seq[best_idx]
|
|
312 query_seq_best = query_seq[best_idx]
|
|
313 domain_size_best = domain_size[best_idx]
|
|
314 [percent_ident, align_similarity, relat_align_len, relat_interrupt,
|
|
315 db_len_proportion
|
|
316 ] = filter_params(db_seq_best, query_seq_best, domain_size_best)
|
|
317 ann_substring = "|".join(ann_substring.split("|")[1:])
|
|
318 annotation_best = "|".join([db_name_best] + annotation_best.split("|")[1:])
|
|
319 if "DANTE_PART" in seq_id:
|
|
320 part = int(seq_id.split("DANTE_PART")[1].split(":")[0].split("_")[0])
|
|
321 dom_start = dom_start + (part - 1) * step
|
|
322 dom_end = dom_end + (part - 1) * step
|
|
323 best_start = best_start + (part - 1) * step
|
|
324 best_end = best_end + (part - 1) * step
|
|
325 if ann_substring is '':
|
|
326 ann_substring = "NONE(Annotations from different classes)"
|
|
327 if len(unique_annotations) > 1:
|
|
328 unique_annotations = ",".join(["{}[{}bp]".format(
|
|
329 ann, pos) for ann, pos in zip(unique_annotations, ann_pos_counts)])
|
|
330 else:
|
|
331 unique_annotations = unique_annotations[0]
|
|
332 if __name__ == '__main__':
|
|
333 SOURCE = configuration.SOURCE_DANTE
|
|
334 else:
|
|
335 SOURCE = configuration.SOURCE_PROFREP
|
|
336 if "/" in domain_type:
|
|
337 gff.write(
|
|
338 "{}\t{}\t{}\t{}\t{}\t.\t{}\t{}\tName={};Final_Classification=Ambiguous_domain;Region_Hits_Classifications_={}\n".format(
|
|
339 seq_id, SOURCE, configuration.DOMAINS_FEATURE, dom_start,
|
|
340 dom_end, strand, configuration.PHASE, domain_type,
|
|
341 unique_annotations))
|
|
342 else:
|
|
343 gff.write(
|
|
344 "{}\t{}\t{}\t{}\t{}\t{}\t{}\t{}\tName={};Final_Classification={};Region_Hits_Classifications={};Best_Hit={}:{}-{}[{}percent];Best_Hit_DB_Pos={}:{}of{};DB_Seq={};Query_Seq={};Identity={};Similarity={};Relat_Length={};Relat_Interruptions={};Hit_to_DB_Length={}\n".format(
|
|
345 seq_id, SOURCE, configuration.DOMAINS_FEATURE, dom_start,
|
|
346 dom_end, best_score, strand, configuration.PHASE, domain_type,
|
|
347 ann_substring, unique_annotations, annotation_best, best_start,
|
|
348 best_end, length_proportion, db_starts_best, db_ends_best,
|
|
349 domain_size_best, db_seq_best, query_seq_best, percent_ident,
|
|
350 align_similarity, relat_align_len, relat_interrupt,
|
|
351 db_len_proportion))
|
|
352
|
|
353
|
|
354 def filter_params(db, query, protein_len):
|
|
355 ''' Calculate basic statistics of the quality of the alignment '''
|
|
356 score_dict = alignment_scoring()
|
|
357 num_ident = 0
|
|
358 count_interrupt = 0
|
|
359 count_similarity = 0
|
|
360 alignment_len = 0
|
|
361 for i, j in zip(db.upper(), query.upper()):
|
|
362 if i == j and i != "X":
|
|
363 num_ident += 1
|
|
364 if j == "/" or j == "\\" or j == "*":
|
|
365 count_interrupt += 1
|
|
366 if (i.isalpha() or i == "*") and (j.isalpha() or j == "*"):
|
|
367 if int(score_dict["{}{}".format(i, j)]) > 0:
|
|
368 count_similarity += 1
|
|
369 ## gapless alignment length proportional to the domain protein length
|
|
370 relat_align_len = round((len(db) - db.count("-")) / protein_len, 3)
|
|
371 ## proportional identical bases (except of X) to al.length
|
|
372 align_identity = round(num_ident / len(db), 2)
|
|
373 ## proportional count of positive scores from scoring matrix to al. length
|
|
374 align_similarity = round(count_similarity / len(db), 2)
|
|
375 ## number of interruptions per 100 bp
|
|
376 relat_interrupt = round(count_interrupt / math.ceil((len(query) / 100)), 2)
|
|
377 ## Proportion of alignment to the original length of protein domain from database (indels included)
|
|
378 db_len_proportion = round(len(db) / protein_len, 2)
|
|
379 return align_identity, align_similarity, relat_align_len, relat_interrupt, db_len_proportion
|
|
380
|
|
381
|
|
382 def line_generator(tab_pipe, maf_pipe, start):
|
|
383 ''' Yield individual lines of LASTAL stdout for single sequence '''
|
|
384 if hasattr(line_generator, "dom"):
|
|
385 seq_id = line_generator.dom.split("\t")[6]
|
|
386 yield line_generator.dom.encode("utf-8")
|
|
387 del line_generator.dom
|
|
388 line_tab = ""
|
|
389 for line_tab in tab_pipe:
|
|
390 line_tab = line_tab.decode("utf-8")
|
|
391 if not line_tab.startswith('#'):
|
|
392 if start:
|
|
393 if not ('seq_id' in locals() and
|
|
394 seq_id != line_tab.split("\t")[6]):
|
|
395 seq_id = line_tab.split("\t")[6]
|
|
396 start = False
|
|
397 line_maf = [maf_pipe.readline() for line_count in range(4)]
|
|
398 db_seq = line_maf[1].decode("utf-8").rstrip().split(" ")[-1]
|
|
399 alignment_seq = line_maf[2].decode("utf-8").rstrip().split(" ")[-1]
|
|
400 line = "{}\t{}\t{}".format(line_tab, db_seq, alignment_seq)
|
|
401 line_id = line.split("\t")[6]
|
|
402 if seq_id != line_id:
|
|
403 line_generator.dom = line
|
|
404 return
|
|
405 else:
|
|
406 yield line.encode("utf-8")
|
|
407 else:
|
|
408 maf_pipe.readline()
|
|
409 if line_tab == "":
|
|
410 raise RuntimeError
|
|
411 else:
|
|
412 return
|
|
413
|
|
414
|
|
415 def get_version(path, LAST_DB):
|
5
|
416 '''Return version is run from git repository '''
|
|
417 try:
|
|
418 branch = subprocess.check_output("git rev-parse --abbrev-ref HEAD",
|
|
419 shell=True,
|
|
420 cwd=path).decode('ascii').strip()
|
|
421 shorthash = subprocess.check_output("git log --pretty=format:'%h' -n 1 ",
|
|
422 shell=True,
|
|
423 cwd=path).decode('ascii').strip()
|
|
424 revcount = len(subprocess.check_output("git log --oneline",
|
|
425 shell=True,
|
|
426 cwd=path).decode('ascii').split())
|
|
427 version_string = (
|
|
428 "##-----------------------------------------------\n"
|
|
429 "##PIPELINE VERSION : "
|
|
430 "{branch}-rv-{revcount}({shorthash})\n"
|
|
431 "##PROTEIN DATABASE VERSION : {PD}\n"
|
|
432 "##-----------------------------------------------\n").format(
|
|
433 branch=branch,
|
|
434 shorthash=shorthash,
|
|
435 revcount=revcount,
|
|
436 PD=os.path.basename(LAST_DB))
|
|
437 except:
|
|
438 version_string = (
|
|
439 "##-----------------------------------------------\n"
|
|
440 "##PROTEIN DATABASE VERSION : {PD}\n"
|
|
441 "##-----------------------------------------------\n").format(
|
|
442 PD=os.path.basename(LAST_DB)
|
|
443 )
|
|
444
|
0
|
445 return version_string
|
|
446
|
|
447
|
|
448 def write_info(dom_gff_tmp, version_string):
|
|
449 dom_gff_tmp.write("{}\n".format(configuration.HEADER_GFF))
|
|
450 dom_gff_tmp.write(version_string)
|
|
451
|
|
452
|
|
453 def domain_search(QUERY, LAST_DB, CLASSIFICATION, OUTPUT_DOMAIN,
|
|
454 THRESHOLD_SCORE, WIN_DOM, OVERLAP_DOM):
|
|
455 ''' Search for protein domains using our protein database and external tool LAST,
|
|
456 stdout is parsed in real time and hits for a single sequence undergo further processing
|
|
457 - tabular format(TAB) to get info about position, score, orientation
|
|
458 - MAF format to gain alignment and original sequence
|
|
459 '''
|
|
460
|
|
461 step = WIN_DOM - OVERLAP_DOM
|
|
462 [headers, above_win, below_win, lens_above_win, seq_starts, seq_ends
|
|
463 ] = characterize_fasta(QUERY, WIN_DOM)
|
|
464 query_temp = split_fasta(QUERY, WIN_DOM, step, headers, above_win,
|
|
465 below_win, lens_above_win, seq_starts, seq_ends)
|
|
466
|
|
467 ## TAB output contains all the alignment scores, positions, strands...
|
|
468 tab = subprocess.Popen(
|
|
469 "lastal -F15 {} {} -L 10 -m 70 -p BL80 -f TAB".format(LAST_DB,
|
|
470 query_temp),
|
|
471 stdout=subprocess.PIPE,
|
|
472 shell=True)
|
|
473 ## MAF output contains alignment sequences
|
|
474 maf = subprocess.Popen(
|
|
475 "lastal -F15 {} {} -L 10 -m 70 -p BL80 -f MAF".format(LAST_DB,
|
|
476 query_temp),
|
|
477 stdout=subprocess.PIPE,
|
|
478 shell=True)
|
|
479
|
|
480 tab_pipe = tab.stdout
|
|
481 maf_pipe = maf.stdout
|
|
482 maf_pipe.readline()
|
|
483
|
|
484 seq_ids = []
|
|
485 dom_tmp = NamedTemporaryFile(delete=False)
|
|
486 with open(dom_tmp.name, "w") as dom_gff_tmp:
|
|
487 path = os.path.dirname(os.path.realpath(__file__))
|
|
488 version_string = get_version(path, LAST_DB)
|
|
489 write_info(dom_gff_tmp, version_string)
|
|
490 gff = open(dom_tmp.name, "a")
|
|
491 start = True
|
|
492 while True:
|
|
493 try:
|
|
494 with warnings.catch_warnings():
|
|
495 warnings.simplefilter("ignore")
|
|
496 sequence_hits = np.genfromtxt(
|
|
497 line_generator(tab_pipe, maf_pipe, start),
|
|
498 names=
|
|
499 "score, name_db, start_db, al_size_db, strand_db, seq_size_db, name_q, start_q, al_size_q, strand_q, seq_size_q, block1, block2, block3, db_seq, q_seq",
|
|
500 usecols=
|
|
501 "score, name_q, start_q, al_size_q, strand_q, seq_size_q, name_db, db_seq, q_seq, seq_size_db, start_db, al_size_db",
|
|
502 dtype=None,
|
|
503 comments=None)
|
|
504 except RuntimeError:
|
|
505 break
|
|
506 ## if there are no domains found
|
|
507 if sequence_hits.size is 0:
|
|
508 gff.write("##NO DOMAINS")
|
|
509 return [], [], [], []
|
|
510
|
|
511 ############# PARSING LASTAL OUTPUT ############################
|
|
512 sequence_hits = np.atleast_1d(sequence_hits)
|
|
513 score = sequence_hits['score'].astype("int")
|
|
514 seq_id = sequence_hits['name_q'][0].astype("str")
|
|
515 start_hit = sequence_hits['start_q'].astype("int")
|
|
516 end_hit = start_hit + sequence_hits['al_size_q'].astype("int")
|
|
517 strand = sequence_hits['strand_q'].astype("str")
|
|
518 seq_len = sequence_hits['seq_size_q'].astype("int")
|
|
519 domain_db = sequence_hits['name_db'].astype("str")
|
|
520 db_seq = sequence_hits['db_seq'].astype("str")
|
|
521 query_seq = sequence_hits['q_seq'].astype("str")
|
|
522 domain_size = sequence_hits['seq_size_db'].astype("int")
|
|
523 db_start = sequence_hits['start_db'].astype("int") + 1
|
|
524 db_end = sequence_hits['start_db'].astype("int") + sequence_hits[
|
|
525 'al_size_db'].astype("int")
|
|
526
|
|
527 [reverse_strand_idx, positions_plus, positions_minus
|
|
528 ] = hits_processing(seq_len, start_hit, end_hit, strand)
|
|
529 strand_gff = "+"
|
|
530 [mins_plus, maxs_plus, data_plus, indices_plus
|
|
531 ] = overlapping_regions(positions_plus)
|
|
532 [mins_minus, maxs_minus, data_minus, indices_minus
|
|
533 ] = overlapping_regions(positions_minus)
|
|
534 positions = positions_plus + positions_minus
|
|
535 indices_overal = indices_plus + [x + reverse_strand_idx
|
|
536 for x in indices_minus]
|
|
537 mins = mins_plus + mins_minus
|
|
538 maxs = maxs_plus + maxs_minus
|
|
539 data = data_plus + data_minus
|
|
540 ## process every region (cluster) of overlapping hits sequentially
|
|
541 count_region = 0
|
|
542 for region in indices_overal:
|
|
543 db_names = domain_db[np.array(region)]
|
|
544 db_starts = db_start[np.array(region)]
|
|
545 db_ends = db_end[np.array(region)]
|
|
546 scores = score[np.array(region)]
|
|
547 annotations = domain_annotation(db_names, CLASSIFICATION)
|
|
548 [score_matrix, classes_dict] = score_table(
|
|
549 mins[count_region], maxs[count_region], data[count_region],
|
|
550 annotations, scores, CLASSIFICATION)
|
|
551 ann_per_reg = score_matrix_evaluation(score_matrix, classes_dict,
|
|
552 THRESHOLD_SCORE)
|
|
553 [domain_type, ann_substring, unique_annotations, ann_pos_counts
|
|
554 ] = group_annot_regs(ann_per_reg)
|
|
555 [best_idx, best_idx_reg] = best_score(scores, region)
|
|
556 annotation_best = annotations[best_idx_reg]
|
|
557 db_name_best = db_names[best_idx_reg]
|
|
558 db_starts_best = db_starts[best_idx_reg]
|
|
559 db_ends_best = db_ends[best_idx_reg]
|
|
560 if count_region == len(indices_plus):
|
|
561 strand_gff = "-"
|
|
562 create_gff3(domain_type, ann_substring, unique_annotations,
|
|
563 ann_pos_counts, mins[count_region], maxs[count_region],
|
|
564 step, best_idx, annotation_best, db_name_best,
|
|
565 db_starts_best, db_ends_best, strand_gff, score,
|
|
566 seq_id, db_seq, query_seq, domain_size, positions, gff)
|
|
567 count_region += 1
|
|
568 seq_ids.append(seq_id)
|
|
569 os.unlink(query_temp)
|
|
570 gff.close()
|
|
571 dom_tmp.close()
|
|
572 ## if any sequence from input data was split into windows, merge and adjust the data from individual windows
|
|
573 if any("DANTE_PART" in x for x in seq_ids):
|
|
574 adjust_gff(OUTPUT_DOMAIN, dom_tmp.name, WIN_DOM, OVERLAP_DOM, step)
|
|
575 ## otherwise use the temporary output as the final domains gff
|
|
576 else:
|
|
577 shutil.copy2(dom_tmp.name, OUTPUT_DOMAIN)
|
|
578 os.unlink(dom_tmp.name)
|
|
579
|
|
580
|
|
581 def adjust_gff(OUTPUT_DOMAIN, gff, WIN_DOM, OVERLAP_DOM, step):
|
|
582 ''' Original gff file is adjusted in case of containing cut parts
|
|
583 - for consecutive sequences overlap is divided to half with first half
|
|
584 of records(domains) belonging to the first sequence and second to the following one.
|
|
585 Duplicate domains going through the middle of the overlap are removed.
|
|
586 First and the last part (marked as LAST) of a certain sequence are
|
|
587 handled separately as the are overlapped from one side only '''
|
|
588
|
|
589 seq_id_all = []
|
|
590 class_dict = defaultdict(int)
|
|
591 seen = set()
|
|
592 with open(OUTPUT_DOMAIN, "w") as adjusted_gff:
|
|
593 with open(gff, "r") as primary_gff:
|
|
594 start = True
|
|
595 for line in primary_gff:
|
|
596 if line.startswith("#"):
|
|
597 adjusted_gff.write(line)
|
|
598 else:
|
|
599 split_line = line.split("\t")
|
|
600 classification = split_line[-1].split(";")[1].split("=")[1]
|
|
601 if start:
|
|
602 seq_id_all.append(split_line[0].split("_DANTE_PART")[
|
|
603 0])
|
|
604 start = False
|
|
605 seq_id = split_line[0].split("_DANTE_PART")[0]
|
|
606 if "DANTE_PART" in line:
|
|
607 line_without_id = "\t".join(split_line[1:])
|
|
608 part = int(split_line[0].split("_DANTE_PART")[1].split(
|
|
609 ":")[0].split("_")[0])
|
|
610 if seq_id != seq_id_all[-1]:
|
|
611 seq_id_all.append(seq_id)
|
|
612
|
|
613 ## first part of the sequence
|
|
614 if part == 1:
|
|
615 cut_end = WIN_DOM - OVERLAP_DOM / 2
|
|
616 if int(split_line[3]) <= cut_end <= int(split_line[
|
|
617 4]):
|
|
618 if line_without_id not in seen:
|
|
619 adjusted_gff.write("{}\t{}".format(
|
|
620 seq_id, line_without_id))
|
|
621 class_dict[classification] += 1
|
|
622 seen.add(line_without_id)
|
|
623 elif int(split_line[4]) < cut_end:
|
|
624 adjusted_gff.write("{}\t{}".format(
|
|
625 seq_id, line_without_id))
|
|
626 class_dict[classification] += 1
|
|
627
|
|
628 ## last part of the sequence
|
|
629 elif "LAST" in split_line[0]:
|
|
630 cut_start = OVERLAP_DOM / 2 + (part - 1) * step
|
|
631 if int(split_line[3]) <= cut_start <= int(
|
|
632 split_line[4]):
|
|
633 if line_without_id not in seen:
|
|
634 adjusted_gff.write("{}\t{}".format(
|
|
635 seq_id, line_without_id))
|
|
636 class_dict[classification] += 1
|
|
637 seen.add(line_without_id)
|
|
638 elif int(split_line[3]) > cut_start:
|
|
639 adjusted_gff.write("{}\t{}".format(
|
|
640 seq_id, line_without_id))
|
|
641 class_dict[classification] += 1
|
|
642
|
|
643 ## middle part of the sequence
|
|
644 else:
|
|
645 cut_start = OVERLAP_DOM / 2 + (part - 1) * step
|
|
646 cut_end = WIN_DOM - OVERLAP_DOM / 2 + (part -
|
|
647 1) * step
|
|
648 if int(split_line[3]) <= cut_start <= int(
|
|
649 split_line[4]) or int(split_line[
|
|
650 3]) <= cut_end <= int(split_line[4]):
|
|
651 if line_without_id not in seen:
|
|
652 adjusted_gff.write("{}\t{}".format(
|
|
653 seq_id, line_without_id))
|
|
654 class_dict[classification] += 1
|
|
655 seen.add(line_without_id)
|
|
656 elif int(split_line[3]) > cut_start and int(
|
|
657 split_line[4]) < cut_end:
|
|
658 adjusted_gff.write("{}\t{}".format(
|
|
659 seq_id, line_without_id))
|
|
660 class_dict[classification] += 1
|
|
661 ## not divived
|
|
662 else:
|
|
663 if seq_id != seq_id_all[-1]:
|
|
664 seq_id_all.append(seq_id)
|
|
665 adjusted_gff.write(line)
|
|
666 class_dict[classification] += 1
|
|
667
|
|
668
|
|
669 def main(args):
|
|
670
|
|
671 t = time.time()
|
|
672
|
|
673 QUERY = args.query
|
|
674 LAST_DB = args.protein_database
|
|
675 CLASSIFICATION = args.classification
|
|
676 OUTPUT_DOMAIN = args.domain_gff
|
|
677 NEW_LDB = args.new_ldb
|
|
678 OUTPUT_DIR = args.output_dir
|
|
679 THRESHOLD_SCORE = args.threshold_score
|
|
680 WIN_DOM = args.win_dom
|
|
681 OVERLAP_DOM = args.overlap_dom
|
|
682
|
|
683 if OUTPUT_DOMAIN is None:
|
|
684 OUTPUT_DOMAIN = configuration.DOMAINS_GFF
|
|
685 if os.path.isdir(LAST_DB):
|
|
686 LAST_DB = os.path.join(LAST_DB, configuration.LAST_DB_FILE)
|
|
687 if os.path.isdir(CLASSIFICATION):
|
|
688 CLASSIFICATION = os.path.join(CLASSIFICATION, configuration.CLASS_FILE)
|
|
689
|
|
690 if NEW_LDB:
|
|
691 subprocess.call("lastdb -p -cR01 {} {}".format(LAST_DB, LAST_DB),
|
|
692 shell=True)
|
|
693
|
|
694 if OUTPUT_DIR and not os.path.exists(OUTPUT_DIR):
|
|
695 os.makedirs(OUTPUT_DIR)
|
|
696
|
|
697 if not os.path.isabs(OUTPUT_DOMAIN):
|
|
698 if OUTPUT_DIR is None:
|
|
699 OUTPUT_DIR = configuration.TMP
|
|
700 if not os.path.exists(OUTPUT_DIR):
|
|
701 os.makedirs(OUTPUT_DIR)
|
|
702 OUTPUT_DOMAIN = os.path.join(OUTPUT_DIR,
|
|
703 os.path.basename(OUTPUT_DOMAIN))
|
|
704 domain_search(QUERY, LAST_DB, CLASSIFICATION, OUTPUT_DOMAIN,
|
|
705 THRESHOLD_SCORE, WIN_DOM, OVERLAP_DOM)
|
|
706
|
|
707 print("ELAPSED_TIME_DOMAINS = {} s".format(time.time() - t))
|
|
708
|
|
709
|
|
710 if __name__ == "__main__":
|
|
711 import argparse
|
|
712 from argparse import RawDescriptionHelpFormatter
|
|
713
|
|
714 class CustomFormatter(argparse.ArgumentDefaultsHelpFormatter,
|
|
715 argparse.RawDescriptionHelpFormatter):
|
|
716 pass
|
|
717
|
|
718 parser = argparse.ArgumentParser(
|
|
719 description=
|
|
720 '''Script performs similarity search on given DNA sequence(s) in (multi)fasta against our protein domains database of all Transposable element for certain group of organisms (Viridiplantae or Metazoans). Domains are subsequently annotated and classified - in case certain domain has multiple annotations assigned, classifation is derived from the common classification level of all of them. Domains search is accomplished engaging LASTAL alignment tool.
|
|
721
|
|
722 DEPENDENCIES:
|
|
723 - python 3.4 or higher with packages:
|
|
724 -numpy
|
|
725 - lastal 744 or higher [http://last.cbrc.jp/]
|
|
726 - configuration.py module
|
|
727
|
|
728 EXAMPLE OF USAGE:
|
|
729
|
|
730 ./protein_domains_pd.py -q PATH_TO_INPUT_SEQ -pdb PATH_TO_PROTEIN_DB -cs PATH_TO_CLASSIFICATION_FILE
|
|
731
|
|
732 When running for the first time with a new database use -nld option allowing lastal to create indexed database files:
|
|
733
|
|
734 -nld True
|
|
735
|
|
736 ''',
|
|
737 epilog="""""",
|
|
738 formatter_class=CustomFormatter)
|
|
739
|
|
740 requiredNamed = parser.add_argument_group('required named arguments')
|
|
741 requiredNamed.add_argument(
|
|
742 "-q",
|
|
743 "--query",
|
|
744 type=str,
|
|
745 required=True,
|
|
746 help=
|
|
747 'input DNA sequence to search for protein domains in a fasta format. Multifasta format allowed.')
|
|
748 requiredNamed.add_argument('-pdb',
|
|
749 "--protein_database",
|
|
750 type=str,
|
|
751 required=True,
|
|
752 help='protein domains database file')
|
|
753 requiredNamed.add_argument('-cs',
|
|
754 '--classification',
|
|
755 type=str,
|
|
756 required=True,
|
|
757 help='protein domains classification file')
|
|
758 parser.add_argument("-oug",
|
|
759 "--domain_gff",
|
|
760 type=str,
|
|
761 help="output domains gff format")
|
|
762 parser.add_argument(
|
|
763 "-nld",
|
|
764 "--new_ldb",
|
|
765 type=str,
|
|
766 default=False,
|
|
767 help=
|
|
768 "create indexed database files for lastal in case of working with new protein db")
|
|
769 parser.add_argument(
|
|
770 "-dir",
|
|
771 "--output_dir",
|
|
772 type=str,
|
|
773 help="specify if you want to change the output directory")
|
|
774 parser.add_argument(
|
|
775 "-thsc",
|
|
776 "--threshold_score",
|
|
777 type=int,
|
|
778 default=80,
|
|
779 help=
|
|
780 "percentage of the best score in the cluster to be tolerated when assigning annotations per base")
|
|
781 parser.add_argument(
|
|
782 "-wd",
|
|
783 "--win_dom",
|
|
784 type=int,
|
|
785 default=10000000,
|
|
786 help="window to process large input sequences sequentially")
|
|
787 parser.add_argument("-od",
|
|
788 "--overlap_dom",
|
|
789 type=int,
|
|
790 default=10000,
|
|
791 help="overlap of sequences in two consecutive windows")
|
|
792
|
|
793 args = parser.parse_args()
|
|
794 main(args)
|