-
Notifications
You must be signed in to change notification settings - Fork 8
/
cobra.py
1887 lines (1678 loc) · 85.9 KB
/
cobra.py
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
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
#!/usr/bin/env python
# Author: LinXing Chen, UC Berkeley
# cobra v1.2.2
# Contig Overlap Based Re-Assembly
# Modification date: Sep 3, 2023
import os
from Bio import SeqIO
from Bio.Seq import reverse_complement
from collections import defaultdict
from Bio.SeqUtils import GC
import argparse
import pysam
from time import strftime
import itertools
parser = argparse.ArgumentParser(description="This script is used to get higher quality (including circular) virus genomes "
"by joining assembled contigs based on their end overlaps.")
requiredNamed = parser.add_argument_group('required named arguments')
requiredNamed.add_argument("-q", "--query", type=str, help="the query contigs file (fasta format), or the query name "
"list (text file, one column).", required=True)
requiredNamed.add_argument("-f", "--fasta", type=str, help="the whole set of assembled contigs (fasta format).", required=True)
requiredNamed.add_argument("-a", "--assembler", type=str, choices=["idba", "megahit", "metaspades"],
help="de novo assembler used, COBRA not tested for others.", required=True)
requiredNamed.add_argument("-mink", "--mink", type=int, help="the min kmer size used in de novo assembly.", required=True)
requiredNamed.add_argument("-maxk", "--maxk", type=int, help="the max kmer size used in de novo assembly.", required=True)
requiredNamed.add_argument("-m", "--mapping", type=str, help="the reads mapping file in sam or bam format.", required=True)
requiredNamed.add_argument("-c", "--coverage", type=str, help="the contig coverage file (two columns divided by tab).",
required=True)
parser.add_argument("-lm", "--linkage_mismatch", type=int, default=2, help="the max read mapping mismatches for "
"determining if two contigs are spanned by "
"paired reads. [2]")
parser.add_argument("-o", "--output", type=str, help="the name of output folder (default = '<query>_COBRA').")
parser.add_argument("-t", "--threads", type=int, default=16, help="the number of threads for blastn. [16]")
parser.add_argument("-v", "--version", action='version', version='cobra v1.2.2')
args = parser.parse_args()
##
#
global one_path_end, link_pair, cov, parsed_linkage, self_circular, two_paths_end, contig2join, contig_checked, \
contig2join_reason, path_circular, path_circular_end, order_all, added_to_contig, header2seq, \
self_circular_non_expected_overlap, all_joined_query, extended_circular_query, extended_partial_query, \
header2len, is_subset_of
#
cov = {} # the coverage of contigs
header2seq = {}
header2len = {}
link_pair = {} # used to save all overlaps between ends
parsed_linkage = set() # Initialize an empty set to store the parsed linkage information
one_path_end = [] # the end of contigs with one potential join
two_paths_end = [] # the end of contigs with two potential joins
self_circular = set()
self_circular_non_expected_overlap = {}
contig2join = {}
contig_checked = {}
contig2join_reason = {}
path_circular_end = set()
path_circular = set()
is_subset_of = {}
extended_circular_query = set()
extended_partial_query = set()
order_all = {}
added_to_contig = {}
all_joined_query = set()
##
# define functions for analyses
def log_info(step, description, line_feed, log_file):
localtime = ' [20' + strftime("%x").rsplit('/',1)[1] + '/' + strftime("%x").rsplit('/',1)[0] + ' ' + strftime("%X") + '] '
print(step + localtime + description, end=line_feed, file=log_file, flush=True)
def determine_file_format(filename):
"""
Determine the format of the input file.
Returns either "fasta" or "txt" based on the content.
"""
with open(filename, 'r') as f:
first_line = f.readline().strip()
if first_line.startswith(">"):
return "fasta"
else:
return "txt"
def contig_name(end_name):
"""
get contig name from end name
"""
if end_name.endswith('_L') or end_name.endswith('_R') or end_name.endswith('_Lrc') or end_name.endswith('_Rrc'):
return end_name.rsplit('_', 1)[0]
else:
return end_name
def get_target(item):
"""
to get the target for next run of joining
"""
suffix = item.rsplit('_', 1)[1]
base_name = item.rsplit('_', 1)[0]
if suffix == 'Lrc':
return base_name + '_Rrc'
elif suffix == 'Rrc':
return base_name + '_Lrc'
elif suffix == 'R':
return base_name + '_L'
else:
return base_name + '_R'
def are_equal_paths(two_paths):
"""
evaluate if two paths are equal, two_paths is a list including two ends, e.g., a_L, b_R
"""
if two_paths[0].rsplit('_', 1)[1].startswith('L') and two_paths[1].rsplit('_', 1)[1].startswith('R'):
return contig_name(two_paths[0]) + '_R' in one_path_end \
and contig_name(two_paths[1]) + '_L' in one_path_end \
and contig_name(link_pair[contig_name(two_paths[0]) + '_R'][0]) == \
contig_name(link_pair[contig_name(two_paths[1]) + '_L'][0])
elif two_paths[0].rsplit('_', 1)[1].startswith('R') and two_paths[1].rsplit('_', 1)[1].startswith('L'):
return contig_name(two_paths[0]) + '_L' in one_path_end \
and contig_name(two_paths[1]) + '_R' in one_path_end \
and contig_name(link_pair[contig_name(two_paths[0]) + '_L'][0]) == \
contig_name(link_pair[contig_name(two_paths[1]) + '_R'][0])
elif two_paths[0].rsplit('_', 1)[1].startswith('R') and two_paths[1].rsplit('_', 1)[1].startswith('R'):
return contig_name(two_paths[0]) + '_L' in one_path_end \
and contig_name(two_paths[1]) + '_L' in one_path_end \
and contig_name(link_pair[contig_name(two_paths[0]) + '_L'][0]) == \
contig_name(link_pair[contig_name(two_paths[1]) + '_L'][0])
elif two_paths[0].rsplit('_', 1)[1].startswith('L') and two_paths[1].rsplit('_', 1)[1].startswith('L'):
return contig_name(two_paths[0]) + '_R' in one_path_end \
and contig_name(two_paths[1]) + '_R' in one_path_end \
and contig_name(link_pair[contig_name(two_paths[0]) + '_R'][0]) == \
contig_name(link_pair[contig_name(two_paths[1]) + '_R'][0])
else:
return False
def the_dominant_one(two_paths):
"""
get the dominant path from two equal paths
"""
if cov[contig_name(two_paths[0])] >= cov[contig_name(two_paths[1])]:
return two_paths[0]
else:
return two_paths[1]
def could_circulate(point, contig, direction):
"""
check if the path is circular with the current contig included
"""
if len(link_pair[point]) == 2: # point is the same as "target" in join_walker
if direction == 'L':
if contig_name(link_pair[point][0]) == contig:
if cov[contig_name(point)] < 1.5 * cov[contig]:
# 2 times is for repeat, but it is too risky, use 1.5 instead (same below)
return link_pair[point][0].endswith('_R') and (contig + '_R', point) in parsed_linkage
else:
return False
elif contig_name(link_pair[point][1]) == contig:
if cov[contig_name(point)] < 1.5 * cov[contig]:
return link_pair[point][1].endswith('_R') and (contig + '_R', point) in parsed_linkage
else:
return False
else:
return False
elif direction == 'R':
if contig_name(link_pair[point][0]) == contig:
if cov[contig_name(point)] < 1.5 * cov[contig]:
return link_pair[point][0].endswith('_L') and (contig + '_L', point) in parsed_linkage
else:
return False
elif contig_name(link_pair[point][1]) == contig:
if cov[contig_name(point)] < 1.5 * cov[contig]:
return link_pair[point][1].endswith('_L') and (contig + '_L', point) in parsed_linkage
else:
return False
else:
return False
else:
return False
elif len(link_pair[point]) == 1:
if direction == 'L':
if contig_name(link_pair[point][0]) == contig:
return link_pair[point][0].endswith('_R') and (contig + '_R', point) in parsed_linkage
elif direction == 'R':
if contig_name(link_pair[point][0]) == contig:
return link_pair[point][0].endswith('_L') and (contig + '_L', point) in parsed_linkage
else:
return False
def the_better_one(two_paths, contig):
"""
Calculate the absolute differences in coverage
"""
diff0 = abs(cov[contig] - cov[contig_name(two_paths[0])])
diff1 = abs(cov[contig] - cov[contig_name(two_paths[1])])
# If diff0 is 0, return two_paths[0]
if diff0 == 0:
return two_paths[0]
else:
# If diff1 is not 0, calculate the ratio and compare
if diff1 != 0:
ratio = diff0 / diff1
if ratio > 1:
return two_paths[1]
elif ratio <= 1:
return two_paths[0]
# If diff1 is 0, return two_paths[1]
else:
return two_paths[1]
def other_end_is_extendable(end, contig):
"""
check if the other end is extendable
"""
if contig_name(end) not in self_circular:
if 0.5 * cov[contig] <= cov[contig_name(end)] <= 2 * cov[contig]:
if get_target(end) in link_pair.keys():
return len(link_pair[get_target(end)]) > 0
else:
return False
else:
return False
else:
return False
def is_ok_to_add(end, contig):
"""
if the other end of a potential joining contig cannot be extended,
add only when it has a very similar coverage to that of the query contig
"""
if contig_name(end) not in self_circular:
return 0.9 * cov[contig] <= cov[contig_name(end)] <= 1.11 * cov[contig]
def not_checked(end_list, checked):
"""
to see if a potential contig has been checked for adding or not
"""
total_checked = 0
for end in end_list:
if contig_name(end) in checked:
total_checked += 1
return total_checked == 0
def detect_self_circular(contig):
"""
to determine the input sequences if they are self_circular genomes or not
"""
end = contig + '_L'
if end in one_path_end:
if link_pair[end][0] == contig + '_R':
# the other end could be joined with the current working end
self_circular.add(contig)
else:
pass
elif end in two_paths_end:
if link_pair[end][0].rsplit('_', 1)[0] != link_pair[end][1].rsplit('_', 1)[0]:
if contig + '_R' in link_pair[end]:
self_circular.add(contig)
else:
pass
else:
pass
else:
pass
def join_walker(contig, direction):
"""
get potential joins for a given query
"""
global contig2join, contig_checked, contig2join_reason, path_circular, path_circular_end
end = contig + '_' + direction
a = len(contig2join[end])
if a == 0:
contig_checked[end].append(contig)
if end in one_path_end:
if contig_name(link_pair[end][0]) != contig:
if other_end_is_extendable(link_pair[end][0], contig) and (end, link_pair[end][0]) in parsed_linkage \
and contig_name(link_pair[end][0]) not in self_circular:
contig2join[end].append(link_pair[end][0])
contig_checked[end].append(contig_name(link_pair[end][0]))
contig2join_reason[contig][contig_name(link_pair[end][0])] = 'other_end_is_extendable'
elif is_ok_to_add(link_pair[end][0], contig) and (end, link_pair[end][0]) in parsed_linkage \
and contig_name(link_pair[end][0]) not in self_circular:
contig2join[end].append(link_pair[end][0])
contig_checked[end].append(contig_name(link_pair[end][0]))
contig2join_reason[contig][contig_name(link_pair[end][0])] = 'is_ok_to_add'
else:
pass
else:
pass
elif end in two_paths_end:
if link_pair[end][0].rsplit('_', 1)[0] != link_pair[end][1].rsplit('_', 1)[0]:
if are_equal_paths(link_pair[end]):
if contig_name(the_dominant_one(link_pair[end])) not in self_circular:
contig2join[end].append(the_dominant_one(link_pair[end]))
contig_checked[end].append(contig_name((link_pair[end][0])))
contig_checked[end].append(contig_name((link_pair[end][1])))
contig2join_reason[contig][contig_name(the_dominant_one(link_pair[end]))] = 'are_equal_paths'
else:
pass
else:
if cov[contig_name(link_pair[end][0])] + cov[contig_name(link_pair[end][1])] >= cov[contig] * 0.5:
# 0.5 is ok, too big will get much fewer "Extended circular" ones.
if contig_name(the_better_one(link_pair[end], contig)) not in self_circular:
contig2join[end].append(the_better_one(link_pair[end], contig))
contig_checked[end].append(contig_name((link_pair[end][0])))
contig_checked[end].append(contig_name((link_pair[end][1])))
contig2join_reason[contig][contig_name(the_better_one(link_pair[end], contig))] = 'the_better_one'
else:
pass
else:
pass
else:
pass
else:
pass
else:
target = get_target(contig2join[end][-1])
if target in one_path_end:
if not_checked(link_pair[target], contig_checked[end]):
if other_end_is_extendable(link_pair[target][0], contig) and (target, link_pair[target][0]) in parsed_linkage \
and contig_name(link_pair[target][0]) not in self_circular:
contig2join[end].append(link_pair[target][0])
contig_checked[end].append(contig_name(link_pair[target][0]))
contig2join_reason[contig][contig_name(link_pair[target][0])] = 'other_end_is_extendable'
elif is_ok_to_add(link_pair[target][0], contig) and (target, link_pair[target][0]) in parsed_linkage \
and contig_name(link_pair[target][0]) not in self_circular:
contig2join[end].append(link_pair[target][0])
contig_checked[end].append(contig_name(link_pair[target][0]))
contig2join_reason[contig][contig_name(link_pair[target][0])] = 'is_ok_to_add'
else:
pass
else:
if could_circulate(target, contig, direction):
path_circular.add(contig)
path_circular_end.add(contig + '_' + direction)
contig2join_reason[contig][contig_name(target)] = 'could_circulate'
elif target in two_paths_end:
if link_pair[target][0].rsplit('_', 1)[0] != link_pair[target][1].rsplit('_', 1)[0]:
if cov[contig_name(target)] >= 1.9 * cov[contig]:
pass
else:
if not_checked(link_pair[target], contig_checked[end]):
if are_equal_paths(link_pair[target]):
if contig_name(the_dominant_one(link_pair[target])) not in self_circular:
contig2join[end].append(the_dominant_one(link_pair[target]))
contig_checked[end].append(contig_name(link_pair[target][0]))
contig_checked[end].append(contig_name(link_pair[target][1]))
contig2join_reason[contig][contig_name(the_dominant_one(link_pair[target]))] = 'are_equal_paths'
else:
pass
else:
if cov[contig_name(link_pair[target][0])] + cov[contig_name(link_pair[target][1])] >= cov[contig] * 0.5:
# 0.5 is ok, too big will get much fewer "Extended circular" ones.
if contig_name(the_better_one(link_pair[target], contig)) not in self_circular:
contig2join[end].append(the_better_one(link_pair[target], contig))
contig_checked[end].append(contig_name((link_pair[target][0])))
contig_checked[end].append(contig_name((link_pair[target][1])))
contig2join_reason[contig][contig_name(the_better_one(link_pair[target], contig))] = 'the_better_one'
else:
pass
else:
pass
else:
if could_circulate(target, contig, direction):
path_circular.add(contig)
path_circular_end.add(contig + '_' + direction)
contig2join_reason[contig][contig_name(target)] = 'could_circulate'
else:
pass
else:
pass
else:
pass
return a < len(contig2join[end])
def join_seqs(contig):
"""
get the join order of the sequences in a give path
"""
global order_all, added_to_contig
order_all[contig] = []
left = contig + '_L'
right = contig + '_R'
added_to_contig[contig] = []
if contig not in path_circular:
if left in contig2join.keys() and right in contig2join.keys():
for item in contig2join[left][::-1]: # the order of contigs added into the left direction should be reversed
order_all[contig].append(item)
added_to_contig[contig].append(contig_name(item))
order_all[contig].append(contig) # the query contig itself should be added to the path as well
for item in contig2join[right]:
if contig_name(item) not in added_to_contig[contig]:
order_all[contig].append(item)
elif left in contig2join.keys() and right not in contig2join.keys():
for item in contig2join[left][::-1]: # the order of contigs added into the left direction should be reversed
order_all[contig].append(item)
order_all[contig].append(contig) # the query contig itself should be added to the path as well
else:
order_all[contig].append(contig) # the query contig itself should be added to the path as well
for item in contig2join[right]:
order_all[contig].append(item)
else:
if left in path_circular_end and right in path_circular_end:
for item in contig2join[left][::-1]: # the order of contigs added into the left direction should be reversed
order_all[contig].append(item)
order_all[contig].append(contig) # the query contig itself should be added to the path as well
elif left in path_circular_end and right not in path_circular_end:
for item in contig2join[left][::-1]: # the order of contigs added into the left direction should be reversed
order_all[contig].append(item)
order_all[contig].append(contig) # the query contig itself should be added to the path as well
elif right in path_circular_end and left not in path_circular_end:
order_all[contig].append(contig) # the query contig itself should be added to the path as well
for item in contig2join[right]:
order_all[contig].append(item)
def retrieve(contig):
"""
to retrieve and save all the contigs in the joining path of a query
"""
if contig in order_all.keys() and len(order_all[contig]) > 0:
out = open('COBRA_retrieved_for_joining/{0}_retrieved.fa'.format(contig), 'w')
added = []
for item in order_all[contig]:
if contig_name(item) not in added:
out.write('>' + contig_name(item) + '\n')
out.write(header2seq[contig_name(item)] + '\n')
added.append(contig_name(item))
else:
pass
out.close()
else:
pass
def count_seq(fasta_file):
"""
calculate the number of seqs in a fasta file
"""
seq_num = 0
a = open(fasta_file, 'r')
for line in a.readlines():
if line.startswith('>'):
seq_num += 1
a.close()
return seq_num
def count_len(fasta_file):
"""
calculate the length of sequences in a fasta file
"""
seq_len = 0
a = open(fasta_file, 'r')
for line in a.readlines():
if not line.startswith('>'):
seq_len += len(line.strip())
a.close()
return seq_len
def summary_fasta(fasta_file):
"""
summary basic information of a fasta file
"""
summary_file = open('{0}.summary.txt'.format(fasta_file), 'w')
if 'self_circular' in fasta_file:
summary_file_headers = ['SeqID', 'Length', 'Coverage', 'GC', 'Ns', 'DTR_length', '\n']
summary_file.write('\t'.join(summary_file_headers[:]))
summary_file.flush()
else:
summary_file_headers = ['SeqID', 'Length', 'Coverage', 'GC', 'Ns', '\n']
summary_file.write('\t'.join(summary_file_headers[:]))
summary_file.flush()
with open('{0}'.format(fasta_file), 'r') as f:
for record in SeqIO.parse(f, "fasta"):
header = str(record.id).strip()
seq = str(record.seq)
ns = seq.count('N')
if header.split('_self')[0] in self_circular:
length = args.maxk - 1 if args.assembler == "idba" else args.maxk
sequence_stats = [header, str(len(seq)), str(cov[header.split('_self')[0]]), str(round(GC(seq), 3)), str(ns), str(length), '\n']
summary_file.write('\t'.join(sequence_stats[:]))
elif header.split('_self')[0] in self_circular_non_expected_overlap.keys():
sequence_stats = [header, str(len(seq)), str(cov[header.split('_self')[0]]), str(round(GC(seq), 3)), str(ns),
str(self_circular_non_expected_overlap[header.split('_self')[0]]), '\n']
summary_file.write('\t'.join(sequence_stats[:]))
else:
sequence_stats = [header, str(len(seq)), str(cov[header.split('_extended')[0]]), str(round(GC(seq), 3)), str(ns), '\n']
summary_file.write('\t'.join(sequence_stats[:]))
f.close()
def get_joined_seqs(fasta_file):
joined_seqs = []
a = open(fasta_file, 'r')
for line in a.readlines():
if line.startswith('>'):
joined_seqs.append(line.strip().split(' ')[0][1:])
else:
pass
for item in joined_seqs:
all_joined_query.add(item)
return ','.join(joined_seqs[:])
def summarize(contig):
"""
summary the retrieved contigs and joined information of each query
"""
if contig not in is_subset_of.keys():
b = count_seq('COBRA_retrieved_for_joining/{0}_retrieved.fa'.format(contig)) # number of retrieved contigs
c = count_len('COBRA_retrieved_for_joining/{0}_retrieved.fa'.format(contig)) # total length of retrieved contigs
d = count_len('COBRA_retrieved_for_joining/{0}_retrieved_joined.fa'.format(contig)) # total length after joining
e = get_joined_seqs('COBRA_retrieved_for_joining/{0}_retrieved.fa'.format(contig))
if contig in extended_circular_query:
return '\t'.join([str(b), e, str(c), str(d), str(d - header2len[contig]), 'Extended_circular'])
elif contig in extended_partial_query:
return '\t'.join([str(b), e, str(c), str(d), str(d - header2len[contig]), 'Extended_partial'])
else:
if is_subset_of[contig] in is_subset_of.keys():
item = is_subset_of[is_subset_of[contig]]
b = count_seq('COBRA_retrieved_for_joining/{0}_retrieved.fa'.format(item)) # number of retrieved contigs
c = count_len('COBRA_retrieved_for_joining/{0}_retrieved.fa'.format(item)) # total length of retrieved contigs
d = count_len('COBRA_retrieved_for_joining/{0}_retrieved_joined.fa'.format(item)) # total length after joining
e = get_joined_seqs('COBRA_retrieved_for_joining/{0}_retrieved.fa'.format(item))
if contig in extended_circular_query:
return '\t'.join([str(b), e, str(c), str(d), str(d - header2len[contig]), 'Extended_circular'])
elif contig in extended_partial_query:
return '\t'.join([str(b), e, str(c), str(d), str(d - header2len[contig]), 'Extended_partial'])
else:
item = is_subset_of[contig]
b = count_seq('COBRA_retrieved_for_joining/{0}_retrieved.fa'.format(item)) # number of retrieved contigs
c = count_len('COBRA_retrieved_for_joining/{0}_retrieved.fa'.format(item)) # total length of retrieved contigs
d = count_len('COBRA_retrieved_for_joining/{0}_retrieved_joined.fa'.format(item)) # total length after joining
e = get_joined_seqs('COBRA_retrieved_for_joining/{0}_retrieved.fa'.format(item))
if contig in extended_circular_query:
return '\t'.join([str(b), e, str(c), str(d), str(d - header2len[contig]), 'Extended_circular'])
elif contig in extended_partial_query:
return '\t'.join([str(b), e, str(c), str(d), str(d - header2len[contig]), 'Extended_partial'])
def get_direction(item):
"""
get the direction in a joining path
"""
if item.endswith('rc'):
return 'reverse'
else:
return 'forward'
def total_length(contig_list):
"""
get the total length of all sequences in a joining path before overlap removing
"""
total = 0
for item in contig_list:
total += header2len[contig_name(item)]
return total
def main():
##
# get information from the input files and parameters and save information
# get the name of the query fasta file
if '/' in args.query:
query_name = '{0}'.format(args.query).rsplit('/', 1)[1]
else:
query_name = '{0}'.format(args.query)
# get the name of the whole contigs fasta file
if '/' in args.fasta:
fasta_name = '{0}'.format(args.fasta).rsplit('/', 1)[1]
else:
fasta_name = '{0}'.format(args.fasta)
# folder of output
if not args.output:
working_dir = '{0}_COBRA'.format(query_name)
else:
working_dir = '{0}'.format(args.output)
# checking if output folder exists
if os.path.exists('{0}'.format(working_dir)):
print('Output folder <{0}> exists, please check.'.format(working_dir))
exit()
else:
os.mkdir('{0}'.format(working_dir))
# determine the length of overlap based on assembler and the largest kmer size
if args.assembler == "idba":
length = args.maxk - 1
else:
length = args.maxk
# write input files information to log file
log = open('{0}/log'.format(working_dir), 'w') # log file
log.write('1. INPUT INFORMATION' + '\n')
log.flush()
if args.assembler == 'idba':
parameters = ['# Assembler: IDBA_UD',
'# Min-kmer: ' + str(args.mink).strip(),
'# Max-kmer: ' + str(args.maxk).strip(),
'# Overlap length: ' + str(length) + ' bp',
'# Read mapping max mismatches for contig linkage: ' + str(args.linkage_mismatch),
'# Query contigs: ' + os.path.abspath(args.query),
'# Whole contig set: ' + os.path.abspath(args.fasta),
'# Mapping file: ' + os.path.abspath(args.mapping),
'# Coverage file: ' + os.path.abspath(args.coverage),
'# Output folder: ' + os.path.abspath(working_dir), '\n']
elif args.assembler == 'metaspades':
parameters = ['# Assembler: metaSPAdes',
'# Min-kmer: ' + str(args.mink).strip(),
'# Max-kmer: ' + str(args.maxk).strip(),
'# Overlap length: ' + str(length) + ' bp',
'# Read mapping max mismatches for contig linkage: ' + str(args.linkage_mismatch),
'# Query contigs: ' + os.path.abspath(args.query),
'# Whole contig set: ' + os.path.abspath(args.fasta),
'# Mapping file: ' + os.path.abspath(args.mapping),
'# Coverage file: ' + os.path.abspath(args.coverage),
'# Output folder: ' + os.path.abspath(working_dir), '\n']
else:
parameters = ['# Assembler: MEGAHIT',
'# Min-kmer: ' + str(args.mink).strip(),
'# Max-kmer: ' + str(args.maxk).strip(),
'# Overlap length: ' + str(length) + ' bp',
'# Read mapping max mismatches for contig linkage: ' + str(args.linkage_mismatch),
'# Query contigs: ' + os.path.abspath(args.query),
'# Whole contig set: ' + os.path.abspath(args.fasta),
'# Mapping file: ' + os.path.abspath(args.mapping),
'# Coverage file: ' + os.path.abspath(args.coverage),
'# Output folder: ' + os.path.abspath(working_dir), '\n']
log.write('\n'.join(parameters[:]))
log.write('2. PROCESSING STEPS' + '\n')
log.flush()
##
# import the whole contigs and save their end sequences
log_info('[01/23]', 'Reading contigs and getting the contig end sequences. ', '', log)
#header2seq = {}
#header2len = {}
gc = {}
L = {}
R = {}
Lrc = {}
Rrc = {}
with open('{0}'.format(args.fasta), 'r') as f:
for record in SeqIO.parse(f, "fasta"):
header = str(record.id).strip()
seq = str(record.seq)
header2seq[header] = seq
gc[header] = str(round(GC(seq), 3))
header2len[header] = len(seq)
L[header + '_L'] = seq[0:length] # the first x bp of left end
Lrc[header + '_Lrc'] = reverse_complement(seq[0:length]) # the reverse sequence of first x bp of left end
R[header + '_R'] = seq[-length:] # the first x bp of right end
Rrc[header + '_Rrc'] = reverse_complement(seq[-length:]) # the reverse sequence of first x bp of right end
log.write('A total of {0} contigs were imported.'.format(len(header2seq.keys())) + '\n')
##
# get potential joins
log_info('[02/23]', 'Getting shared contig ends.', '\n', log)
#link_pair = {} # used to save all overlaps between ends
d_L = defaultdict(set)
d_Lrc = defaultdict(set)
d_R = defaultdict(set)
d_Rrc = defaultdict(set)
for k, v in L.items(): # save header2seq in dictionary with seqs as keys
d_L[v].add(k)
for k, v in Lrc.items():
d_Lrc[v].add(k)
for k, v in R.items():
d_R[v].add(k)
for k, v in Rrc.items():
d_Rrc[v].add(k)
d_L_d_Lrc_shared = set(d_L.keys()).intersection(set(d_Lrc.keys()))
# get the shared seqs between direction pairs (L/Lrc, Lrc/L, L/R, R/L, R/Rrc, Rrc/R, Lrc/Rrc, Rrc/Lrc)
d_L_d_R_shared = set(d_L.keys()).intersection(set(d_R.keys()))
# the d_R_d_L_shared will be included below
d_R_d_Rrc_shared = set(d_R.keys()).intersection(set(d_Rrc.keys()))
d_Rrc_d_Lrc_shared = set(d_Rrc.keys()).intersection(set(d_Lrc.keys()))
##
# get link_pair between ends
for end in d_L_d_Lrc_shared:
for left in d_L[end]: # left is a seq name
for left_rc in d_Lrc[end]: # left_rc is a seq name
if left not in link_pair.keys():
link_pair[left] = [left_rc]
else:
link_pair[left].append(left_rc)
for left_rc in d_Lrc[end]:
for left in d_L[end]:
if left_rc not in link_pair.keys():
link_pair[left_rc] = [left]
else:
link_pair[left_rc].append(left)
for end in d_L_d_R_shared:
for left in d_L[end]:
for right in d_R[end]:
if left not in link_pair.keys():
link_pair[left] = [right]
else:
link_pair[left].append(right)
for right in d_R[end]:
for left in d_L[end]:
if right not in link_pair.keys():
link_pair[right] = [left]
else:
link_pair[right].append(left)
for end in d_R_d_Rrc_shared:
for right in d_R[end]:
for right_rc in d_Rrc[end]:
if right not in link_pair.keys():
link_pair[right] = [right_rc]
else:
link_pair[right].append(right_rc)
for right_rc in d_Rrc[end]:
for right in d_R[end]:
if right_rc not in link_pair.keys():
link_pair[right_rc] = [right]
else:
link_pair[right_rc].append(right)
for end in d_Rrc_d_Lrc_shared:
for right_rc in d_Rrc[end]:
for left_rc in d_Lrc[end]:
if right_rc not in link_pair.keys():
link_pair[right_rc] = [left_rc]
else:
link_pair[right_rc].append(left_rc)
for left_rc in d_Lrc[end]:
for right_rc in d_Rrc[end]:
if left_rc not in link_pair.keys():
link_pair[left_rc] = [right_rc]
else:
link_pair[left_rc].append(right_rc)
##
# save all paired links to a file
log_info('[03/23]', 'Writing contig end joining pairs.', '\n', log)
p = open('{0}/COBRA_end_joining_pairs.txt'.format(working_dir), 'w')
#one_path_end = [] # the end of contigs with one potential join
#two_paths_end = [] # the end of contigs with two potential joins
for item in link_pair.keys():
for point in link_pair[item]:
p.write(item + '\t' + point + '\n') # print link pairs into a file for check if interested
if len(link_pair[item]) == 1: # and len(link_pair[link_pair[item][0]]) > 1:
one_path_end.append(item) # add one joining end to a list, its pair may have one or more joins
elif len(link_pair[item]) == 2 and len(link_pair[link_pair[item][0]]) == 1 and len(
link_pair[link_pair[item][1]]) == 1:
two_paths_end.append(item) # add two joining end to a list, each of its pairs should only have one join
else:
pass
p.close()
##
# read and save the coverage of all contigs
log_info('[04/23]', 'Getting contig coverage information.', '\n', log)
#cov = {}
coverage = open('{0}'.format(args.coverage), 'r')
for line in coverage.readlines():
line = line.strip().split('\t')
cov[line[0]] = round(float(line[1]), 3)
coverage.close()
if len(cov.keys()) < len(header2seq.keys()):
print('Some contigs do not have coverage information. Please check. COBRA exits.')
exit()
else:
pass
##
# open the query file and save the information
log_info('[05/23]', 'Getting query contig list. ', '', log)
query_set = set()
orphan_end_query = set()
non_orphan_end_query = set()
with open('{0}'.format(args.query), 'r') as query_file:
if determine_file_format(args.query) == 'fasta': # if the query file is in fasta format
for record in SeqIO.parse(query_file, 'fasta'):
header = str(record.id).strip()
if header in header2seq.keys(): # some queries may not in the whole assembly, should be rare though.
query_set.add(header)
else:
print('Query {0} is not in your whole contig fasta file, please check!'.format(header), flush=True)
else: # if the query file is in text format
for line in query_file:
header = line.strip().split(' ')[0]
if header in header2seq.keys(): # some queries may not in the whole assembly, should be rare though.
query_set.add(header)
else:
print('Query {0} is not in your whole contig fasta file, please check!'.format(header), flush=True)
# distinguish orphan_end_query and non_orphan_end_query:
for header in query_set:
if header + '_L' not in link_pair.keys() and header + '_R' not in link_pair.keys():
orphan_end_query.add(header)
else:
non_orphan_end_query.add(header)
#
log.write('A total of {0} query contigs were imported.'.format(len(query_set)) + '\n')
log.flush()
##
# get the linkage of contigs based on paired-end reads mapping
log_info('[06/23]', 'Getting contig linkage based on sam/bam. Be patient, this may take long.', '\n', log)
linkage = defaultdict(set) # Initialize a defaultdict to store linked contigs
contig_spanned_by_PE_reads = {} # Initialize a dictionary to store paired-end reads spanning contigs
for contig in orphan_end_query:
contig_spanned_by_PE_reads[contig] = defaultdict(list) # Create a defaultdict(list) for each contig in header2seq
with pysam.AlignmentFile('{0}'.format(args.mapping), 'rb') as map_file:
for line in map_file:
if not line.is_unmapped and line.get_tag("NM") <= args.linkage_mismatch:
# mismatch should not be more than the defined threshold
if line.reference_name != line.next_reference_name:
# Check if the read and its mate map to different contigs
if header2len[line.reference_name] > 1000:
# If the contig length is greater than 1000, determine if the read maps to the left or right end
if line.reference_start <= 500:
linkage[line.query_name].add(line.reference_name + '_L') # left end
else:
linkage[line.query_name].add(line.reference_name + '_R') # right end
else:
# If the contig length is 1000 or less, add both the left and right ends to the linkage
linkage[line.query_name].add(line.reference_name + '_L')
linkage[line.query_name].add(line.reference_name + '_R')
else:
# If the read and its mate map to the same contig, store the read mapped position (start)
if line.reference_name in orphan_end_query:
if line.reference_start <= 500 or header2len[line.reference_name] - line.reference_start <= 500:
contig_spanned_by_PE_reads[line.reference_name][line.query_name].append(
line.reference_start)
else:
pass
else:
pass
else:
pass
#
log_info('[07/23]', 'Parsing the linkage information.', '\n', log)
#parsed_linkage = set() # Initialize an empty set to store the parsed linkage information
for read in linkage.keys():
if len(linkage[read]) >= 2: # Process only reads linked to at least two contigs
# for item in linkage[read]:
for item, item_1 in itertools.combinations(linkage[read], 2):
# Generate unique pairs of linked contigs for the current read using itertools.combinations
if item.rsplit('_', 1)[1] != item_1.rsplit('_', 1)[1]:
# If the contigs have different ends (_L or _R), add the combinations to the parsed_linkage
parsed_linkage.add((item, item_1))
parsed_linkage.add((item_1, item))
parsed_linkage.add((item + 'rc', item_1 + 'rc'))
parsed_linkage.add((item_1 + 'rc', item + 'rc'))
else:
# If the contigs have the same ends, add the combinations with reverse-complement (_rc) to the parsed_linkage
parsed_linkage.add((item, item_1 + 'rc'))
parsed_linkage.add((item_1 + 'rc', item))
parsed_linkage.add((item + 'rc', item_1))
parsed_linkage.add((item_1, item + 'rc'))
else:
pass
linkage = None # remove it as no longer used later
#
parsed_contig_spanned_by_PE_reads = set() # Initialize a set to store the contig spanned by paired-end reads
for contig in orphan_end_query:
for PE in contig_spanned_by_PE_reads[
contig].keys(): # Check if the count is 0 and the contig has exactly two paired-end reads
if len(contig_spanned_by_PE_reads[contig][PE]) == 2:
# Check if the absolute difference between the positions of the two paired-end reads is greater than or equal to
# the length of contig minus 1000 bp
if abs(contig_spanned_by_PE_reads[contig][PE][0] - contig_spanned_by_PE_reads[contig][PE][1]) >= \
header2len[contig] - 1000:
parsed_contig_spanned_by_PE_reads.add(contig)
else:
pass
else:
pass
##
#
log_info('[08/23]', 'Detecting self_circular contigs. ', '\n', log)
#self_circular = set()
for contig in non_orphan_end_query:
detect_self_circular(contig)
debug = open('{0}/debug.txt'.format(working_dir), 'w')
# for debug
print('self_circular', file=debug, flush=True)
print(self_circular, file=debug, flush=True)
# orphan end queries info
# determine potential self_circular contigs from contigs with orphan end
#self_circular_non_expected_overlap = {}
min_over_len = 0
if args.assembler == 'idba':
min_over_len = args.mink - 1
else:
min_over_len = args.mink
for contig in orphan_end_query: # determine if there is DTR for those query with orphan ends, if yes, assign as self_circular as well
if contig in parsed_contig_spanned_by_PE_reads:
sequence = header2seq[contig]
end_part = sequence[-min_over_len:]
if sequence.count(end_part) == 2:
expected_end = sequence.split(end_part)[0] + end_part
if sequence.endswith(expected_end):
self_circular_non_expected_overlap[contig] = len(expected_end)
else:
pass
else:
pass
else:
pass
for contig in self_circular_non_expected_overlap.keys():
orphan_end_query.remove(contig)
# debug
print(self_circular_non_expected_overlap, file=debug, flush=True)
##
# walk the joins
log_info('[09/23]', 'Detecting joins of contigs. ', '', log)
#contig2join = {}
#contig_checked = {}
#contig2join_reason = {}
for contig in query_set:
contig2join[contig + '_L'] = []
contig2join[contig + '_R'] = []
contig_checked[contig + '_L'] = []
contig_checked[contig + '_R'] = []
contig2join_reason[contig] = {}
contig2join_reason[contig][contig] = 'query'