forked from lh3/minimap2
-
Notifications
You must be signed in to change notification settings - Fork 0
/
paftools.js
executable file
·2568 lines (2421 loc) · 81.8 KB
/
paftools.js
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 k8
var paftools_version = '2.17-r982-dirty';
/*****************************
***** Library functions *****
*****************************/
/*******************************
* Command line option parsing *
*******************************/
var getopt = function(args, ostr) {
var oli; // option letter list index
if (typeof(getopt.place) == 'undefined')
getopt.ind = 0, getopt.arg = null, getopt.place = -1;
if (getopt.place == -1) { // update scanning pointer
if (getopt.ind >= args.length || args[getopt.ind].charAt(getopt.place = 0) != '-') {
getopt.place = -1;
return null;
}
if (getopt.place + 1 < args[getopt.ind].length && args[getopt.ind].charAt(++getopt.place) == '-') { // found "--"
++getopt.ind;
getopt.place = -1;
return null;
}
}
var optopt = args[getopt.ind].charAt(getopt.place++); // character checked for validity
if (optopt == ':' || (oli = ostr.indexOf(optopt)) < 0) {
if (optopt == '-') return null; // if the user didn't specify '-' as an option, assume it means null.
if (getopt.place < 0) ++getopt.ind;
return '?';
}
if (oli+1 >= ostr.length || ostr.charAt(++oli) != ':') { // don't need argument
getopt.arg = null;
if (getopt.place < 0 || getopt.place >= args[getopt.ind].length) ++getopt.ind, getopt.place = -1;
} else { // need an argument
if (getopt.place >= 0 && getopt.place < args[getopt.ind].length)
getopt.arg = args[getopt.ind].substr(getopt.place);
else if (args.length <= ++getopt.ind) { // no arg
getopt.place = -1;
if (ostr.length > 0 && ostr.charAt(0) == ':') return ':';
return '?';
} else getopt.arg = args[getopt.ind]; // white space
getopt.place = -1;
++getopt.ind;
}
return optopt;
}
/***********************
* Interval operations *
***********************/
Interval = {};
Interval.sort = function(a)
{
if (typeof a[0] == 'number')
a.sort(function(x, y) { return x - y });
else a.sort(function(x, y) { return x[0] != y[0]? x[0] - y[0] : x[1] - y[1] });
}
Interval.merge = function(a, sorted)
{
if (typeof sorted == 'undefined') sorted = true;
if (!sorted) Interval.sort(a);
var k = 0;
for (var i = 1; i < a.length; ++i) {
if (a[k][1] >= a[i][0])
a[k][1] = a[k][1] > a[i][1]? a[k][1] : a[i][1];
else a[++k] = a[i].slice(0);
}
a.length = k + 1;
}
Interval.index_end = function(a, sorted)
{
if (a.length == 0) return;
if (typeof sorted == 'undefined') sorted = true;
if (!sorted) Interval.sort(a);
a[0].push(0);
var k = 0, k_en = a[0][1];
for (var i = 1; i < a.length; ++i) {
if (k_en <= a[i][0]) {
for (++k; k < i; ++k)
if (a[k][1] > a[i][0])
break;
k_en = a[k][1];
}
a[i].push(k);
}
}
Interval.find_intv = function(a, x)
{
var left = -1, right = a.length;
if (typeof a[0] == 'number') {
while (right - left > 1) {
var mid = left + ((right - left) >> 1);
if (a[mid] > x) right = mid;
else if (a[mid] < x) left = mid;
else return mid;
}
} else {
while (right - left > 1) {
var mid = left + ((right - left) >> 1);
if (a[mid][0] > x) right = mid;
else if (a[mid][0] < x) left = mid;
else return mid;
}
}
return left;
}
Interval.find_ovlp = function(a, st, en)
{
if (a.length == 0 || st >= en) return [];
var l = Interval.find_intv(a, st);
var k = l < 0? 0 : a[l][a[l].length - 1];
var b = [];
for (var i = k; i < a.length; ++i) {
if (a[i][0] >= en) break;
else if (st < a[i][1])
b.push(a[i]);
}
return b;
}
/**********************************
* Reverse and reverse complement *
**********************************/
function fasta_read(fn)
{
var h = {}, gt = '>'.charCodeAt(0);
var file = fn == '-'? new File() : new File(fn);
var buf = new Bytes(), seq = null, name = null, seqlen = [];
while (file.readline(buf) >= 0) {
if (buf[0] == gt) {
if (seq != null && name != null) {
seqlen.push([name, seq.length]);
h[name] = seq;
name = seq = null;
}
var m, line = buf.toString();
if ((m = /^>(\S+)/.exec(line)) != null) {
name = m[1];
seq = new Bytes();
}
} else seq.set(buf);
}
if (seq != null && name != null) {
seqlen.push([name, seq.length]);
h[name] = seq;
}
buf.destroy();
file.close();
return [h, seqlen];
}
function fasta_free(fa)
{
for (var name in fa)
fa[name].destroy();
}
Bytes.prototype.reverse = function()
{
for (var i = 0; i < this.length>>1; ++i) {
var tmp = this[i];
this[i] = this[this.length - i - 1];
this[this.length - i - 1] = tmp;
}
}
// reverse complement a DNA string
Bytes.prototype.revcomp = function()
{
if (Bytes.rctab == null) {
var s1 = 'WSATUGCYRKMBDHVNwsatugcyrkmbdhvn';
var s2 = 'WSTAACGRYMKVHDBNwstaacgrymkvhdbn';
Bytes.rctab = [];
for (var i = 0; i < 256; ++i) Bytes.rctab[i] = 0;
for (var i = 0; i < s1.length; ++i)
Bytes.rctab[s1.charCodeAt(i)] = s2.charCodeAt(i);
}
for (var i = 0; i < this.length>>1; ++i) {
var tmp = this[this.length - i - 1];
this[this.length - i - 1] = Bytes.rctab[this[i]];
this[i] = Bytes.rctab[tmp];
}
if (this.length&1)
this[this.length>>1] = Bytes.rctab[this[this.length>>1]];
}
/********************
***** paftools *****
********************/
/*****************
* Miscellaneous *
*****************/
// liftover
function paf_liftover(args)
{
function read_bed(fn, to_merge)
{
if (fn == null) return null;
if (typeof to_merge == 'undefined') to_merge = true;
var file = fn == '-'? new File() : new File(fn);
var buf = new Bytes();
var bed = {};
while (file.readline(buf) >= 0) {
var t = buf.toString().split("\t");
if (bed[t[0]] == null) bed[t[0]] = [];
bed[t[0]].push([parseInt(t[1]), parseInt(t[2])]);
}
buf.destroy();
file.close();
for (var chr in bed) {
Interval.sort(bed[chr]);
if (to_merge)
Interval.merge(bed[chr], true);
Interval.index_end(bed[chr], true);
}
return bed;
}
var re_cigar = /(\d+)([MID])/g, re_tag = /^(\S\S):([AZif]):(\S+)$/;
var c, to_merge = false, min_mapq = 5, min_len = 50000, max_div = 2.0;
var re = /(\d+)([MID])/g;
while ((c = getopt(args, "mq:l:d:")) != null) {
if (c == 'm') to_merge = true;
else if (c == 'q') min_mapq = parseInt(getopt.arg);
else if (c == 'l') min_len = parseInt(getopt.arg);
else if (c == 'd') max_div = parseFloat(getopt.arg);
}
if (args.length - getopt.ind < 2) {
print("Usage: paftools.js liftover [options] <aln.paf> <query.bed>");
print("Options:");
print(" -q INT min mapping quality [" + min_mapq + "]");
print(" -l INT min alignment length [" + min_len + "]");
print(" -d FLOAT max sequence divergence (>=1 to disable) [1]");
exit(1);
}
var bed = read_bed(args[getopt.ind+1], to_merge);
var file = args[getopt.ind] == '-'? new File() : new File(args[getopt.ind]);
var buf = new Bytes();
while (file.readline(buf) >= 0) {
var t = buf.toString().split("\t");
if (bed[t[0]] == null) continue; // sequence not present in BED; skip
// parse tp and cg tags
var m, tp = null, cg = null;
for (var i = 12; i < t.length; ++i) {
if ((m = re_tag.exec(t[i])) != null) {
if (m[1] == 'tp') tp = m[3];
else if (m[1] == 'cg') cg = m[3];
}
}
if (tp != 'P' && tp != 'I') continue; // only process primary alignments
if (cg == null) throw Error("unable to find the 'cg' tag");
// filter out bad alignments and check overlaps
for (var i = 1; i <= 3; ++i)
t[i] = parseInt(t[i]);
for (var i = 6; i <= 11; ++i)
t[i] = parseInt(t[i]);
if (t[11] < min_mapq || t[10] < min_len) continue;
var regs = Interval.find_ovlp(bed[t[0]], t[2], t[3]);
if (regs.length == 0) continue; // not overlapping any regions in input BED
if (max_div >= 0.0 && max_div < 1.0) {
var n_gaps = 0, n_opens = 0;
while ((m = re_cigar.exec(cg)) != null)
if (m[2] == 'I' || m[2] == 'D')
n_gaps += parseInt(m[1]), ++n_opens;
var n_mm = t[10] - t[9] - n_gaps;
var n_diff2 = n_mm + n_opens;
if (n_diff2 / (n_diff2 + t[9]) > max_div)
continue;
}
// extract start and end positions
var a = [], r = [], strand = t[4];
for (var i = 0; i < regs.length; ++i) {
var s = regs[i][0], e = regs[i][1];
if (strand == '+') {
a.push([s, 0, i, -2]);
a.push([e - 1, 1, i, -2]);
} else {
a.push([t[1] - e, 0, i, -2]);
a.push([t[1] - s - 1, 1, i, -2]);
}
r.push([-2, -2]);
}
a.sort(function(x, y) { return x[0] - y[0] });
// lift start/end positions
var k = 0, x = t[7], y = strand == '+'? t[2] : t[1] - t[3];
while ((m = re_cigar.exec(cg)) != null) { // TODO: be more careful about edge cases
var len = parseInt(m[1]);
if (m[2] == 'D') { // do nothing for D
x += len;
continue;
}
while (k < a.length && a[k][0] < y) ++k; // skip out-of-range positions
for (var i = k; i < a.length; ++i) {
if (y <= a[i][0] && a[i][0] < y + len)
a[i][3] = m[2] == 'M'? x + (a[i][0] - y) : x;
else break;
}
y += len;
if (m[2] == 'M') x += len;
}
if (x != t[8] || (strand == '+' && y != t[3]) || (strand == '-' && y != t[1] - t[2]))
throw Error("CIGAR is inconsistent with mapping coordinates");
// generate result
for (var i = 0; i < a.length; ++i) {
if (a[i][1] == 0) r[a[i][2]][0] = a[i][3];
else r[a[i][2]][1] = a[i][3] + 1; // change to half-close-half-open
}
for (var i = 0; i < r.length; ++i) {
var name = [t[0], regs[i][0], regs[i][1]].join("_");
if (r[i][0] < 0) name += "_t5", r[i][0] = t[7];
if (r[i][1] < 0) name += "_t3", r[i][1] = t[8];
print(t[5], r[i][0], r[i][1], name, 0, strand);
}
}
buf.destroy();
file.close();
}
// variant calling
function paf_call(args)
{
var re_cs = /([:=*+-])(\d+|[A-Za-z]+)/g, re_tag = /\t(\S\S:[AZif]):(\S+)/g;
var c, min_cov_len = 10000, min_var_len = 50000, gap_thres = 50, gap_thres_long = 1000, min_mapq = 5;
var fa_tmp = null, fa, fa_lens, is_vcf = false, sample_name = "sample";
while ((c = getopt(args, "l:L:g:q:B:f:s:")) != null) {
if (c == 'l') min_cov_len = parseInt(getopt.arg);
else if (c == 'L') min_var_len = parseInt(getopt.arg);
else if (c == 'g') gap_thres = parseInt(getopt.arg);
else if (c == 'G') gap_thres_long = parseInt(getopt.arg);
else if (c == 'q') min_mapq = parseInt(getopt.arg);
else if (c == 'f') fa_tmp = fasta_read(getopt.arg, fa_lens);
else if (c == 's') sample_name = getopt.arg;
}
if (fa_tmp != null) fa = fa_tmp[0], fa_lens = fa_tmp[1], is_vcf = true;
if (args.length == getopt.ind) {
print("Usage: sort -k6,6 -k8,8n <with-cs.paf> | paftools.js call [options] -");
print("Options:");
print(" -l INT min alignment length to compute coverage ["+min_cov_len+"]");
print(" -L INT min alignment length to call variants ["+min_var_len+"]");
print(" -q INT min mapping quality ["+min_mapq+"]");
print(" -g INT short/long gap threshold (for statistics only) ["+gap_thres+"]");
print(" -f FILE reference sequences (enabling VCF output) [null]");
print(" -s NAME sample name in VCF header ["+sample_name+"]");
exit(1);
}
var file = args[getopt.ind] == '-'? new File() : new File(args[getopt.ind]);
var buf = new Bytes();
var tot_len = 0, n_sub = [0, 0, 0], n_ins = [0, 0, 0, 0, 0], n_del = [0, 0, 0, 0, 0];
function print_vcf(o, fa)
{
var v = null;
if (o[3] != 1) return; // coverage is one; skip
if (o[5] == '-' && o[6] == '-') return;
if (o[5] != '-' && o[6] != '-') { // snp
v = [o[0], o[1] + 1, '.', o[5].toUpperCase(), o[6].toUpperCase()];
} else if (o[1] > 0) { // shouldn't happen in theory
if (fa[o[0]] == null) throw Error('sequence "' + o[0] + '" is absent from the reference FASTA');
if (o[1] >= fa[o[0]].length) throw Error('position ' + o[1] + ' exceeds the length of sequence "' + o[0] + '"');
var ref = String.fromCharCode(fa[o[0]][o[1]-1]).toUpperCase();
if (o[5] == '-') // insertion
v = [o[0], o[1], '.', ref, ref + o[6].toUpperCase()];
else // deletion
v = [o[0], o[1], '.', ref + o[5].toUpperCase(), ref];
}
v.push(o[4], '.', 'QNAME=' + o[7] + ';QSTART=' + (o[8]+1) + ';QSTRAND=' + (rev? '-' : '+'), 'GT', '1/1');
if (v == null) throw Error("unexpected variant: [" + o.join(",") + "]");
print(v.join("\t"));
}
function count_var(o)
{
if (o[3] > 1) return;
if (o[5] == '-' && o[6] == '-') return;
if (o[5] == '-') { // insertion
var l = o[6].length;
if (l == 1) ++n_ins[0];
else if (l == 2) ++n_ins[1];
else if (l < gap_thres) ++n_ins[2];
else if (l < gap_thres_long) ++n_ins[3];
else ++n_ins[4];
} else if (o[6] == '-') { // deletion
var l = o[5].length;
if (l == 1) ++n_del[0];
else if (l == 2) ++n_del[1];
else if (l < gap_thres) ++n_del[2];
else if (l < gap_thres_long) ++n_del[3];
else ++n_del[4];
} else {
++n_sub[0];
var s = (o[5] + o[6]).toLowerCase();
if (s == 'ag' || s == 'ga' || s == 'ct' || s == 'tc')
++n_sub[1];
else ++n_sub[2];
}
}
if (is_vcf) {
print('##fileformat=VCFv4.1');
for (var i = 0; i < fa_lens.length; ++i)
print('##contig=<ID=' + fa_lens[i][0] + ',length=' + fa_lens[i][1] + '>');
print('##INFO=<ID=QNAME,Number=1,Type=String,Description="Query name">');
print('##INFO=<ID=QSTART,Number=1,Type=Integer,Description="Query start">');
print('##INFO=<ID=QSTRAND,Number=1,Type=String,Description="Query strand">');
print('##FORMAT=<ID=GT,Number=1,Type=String,Description="Genotype">');
print('#CHROM POS ID REF ALT QUAL FILTER INFO FORMAT '+sample_name);
}
var a = [], out = [];
var c1_ctg = null, c1_start = 0, c1_end = 0, c1_counted = false, c1_len = 0;
while (file.readline(buf) >= 0) {
var line = buf.toString();
var m, t = line.split("\t", 12);
if (t.length < 12 || t[5] == '*') continue; // unmapped
for (var i = 6; i <= 11; ++i)
t[i] = parseInt(t[i]);
if (t[10] < min_cov_len || t[11] < min_mapq) continue;
//print(t[0], t[7], t[8], c1_start, c1_end);
for (var i = 1; i <= 3; ++i)
t[i] = parseInt(t[i]);
var ctg = t[5], x = t[7], end = t[8];
var query = t[0], rev = (t[4] == '-'), y = rev? t[3] : t[2];
// collect tags
var cs = null, tp = null, have_s1 = false, have_s2 = false;
while ((m = re_tag.exec(line)) != null) {
if (m[1] == 'cs:Z') cs = m[2];
else if (m[1] == 'tp:A') tp = m[2];
else if (m[1] == 's1:i') have_s1 = true;
else if (m[1] == 's2:i') have_s2 = true;
}
if (have_s1 && !have_s2) continue;
if (tp != null && (tp == 'S' || tp == 'i')) continue;
// compute regions covered by 1 contig
if (ctg != c1_ctg || x >= c1_end) {
if (c1_counted && c1_end > c1_start) {
c1_len += c1_end - c1_start;
if (!is_vcf) print('R', c1_ctg, c1_start, c1_end);
}
c1_ctg = ctg, c1_start = x, c1_end = end;
c1_counted = (t[10] >= min_var_len);
} else if (end > c1_end) { // overlap
if (c1_counted && x > c1_start) {
c1_len += x - c1_start;
if (!is_vcf) print('R', c1_ctg, c1_start, x);
}
c1_start = c1_end, c1_end = end;
c1_counted = (t[10] >= min_var_len);
} else if (end > c1_start) { // contained
if (c1_counted && x > c1_start) {
c1_len += x - c1_start;
if (!is_vcf) print('R', c1_ctg, c1_start, x);
}
c1_start = end;
} // else, the alignment precedes the cov1 region; do nothing
// output variants ahead of this alignment
while (out.length) {
if (out[0][0] != ctg || out[0][2] <= x) {
count_var(out[0]);
if (is_vcf) print_vcf(out[0], fa);
else print('V', out[0].join("\t"));
out.shift();
} else break;
}
// update coverage
for (var i = 0; i < out.length; ++i)
if (out[i][1] >= x && out[i][2] <= end)
++out[i][3];
// drop alignments that don't overlap with the current one
var k = 0;
for (var i = 0; i < a.length; ++i)
if (a[i][0] == ctg && a[i][2] > x)
a[k++] = a[i];
a.length = k;
// core loop
if (t[10] >= min_var_len) {
if (cs == null) continue; // no cs tag
var blen = 0, n_diff = 0;
tot_len += t[10];
while ((m = re_cs.exec(cs)) != null) {
var cov = 1;
if (m[1] == '*' || m[1] == '+' || m[1] == '-')
for (var i = 0; i < a.length; ++i)
if (a[i][2] > x) ++cov;
var qs, qe;
if (m[1] == '=' || m[1] == ':') {
var l = m[1] == '='? m[2].length : parseInt(m[2]);
if (rev) y -= l;
else y += l;
x += l, blen += l;
} else if (m[1] == '*') {
if (rev) qs = y - 1, qe = y, --y;
else qs = y, qe = y + 1, ++y;
var br = m[2].charAt(0), bq = m[2].charAt(1);
if (br != 'n' && bq != 'n') { // don't call a SNP if there is an ambiguous base
out.push([t[5], x, x+1, cov, t[11], br, bq, query, qs, qe, rev? '-' : '+']);
++n_diff;
}
++x, ++blen;
} else if (m[1] == '+') {
var l = m[2].length;
if (rev) qs = y - l, qe = y, y -= l;
else qs = y, qe = y + l, y += l;
out.push([t[5], x, x, cov, t[11], '-', m[2], query, qs, qe, rev? '-' : '+']);
++blen, ++n_diff;
} else if (m[1] == '-') {
var l = m[2].length;
out.push([t[5], x, x + l, cov, t[11], m[2], '-', query, y, y, rev? '-' : '+']);
x += l, ++blen, ++n_diff;
}
}
}
a.push([t[5], t[7], t[8]]);
}
if (c1_counted && c1_end > c1_start) {
c1_len += c1_end - c1_start;
if (!is_vcf) print('R', c1_ctg, c1_start, c1_end);
}
while (out.length) {
count_var(out[0]);
if (is_vcf) print_vcf(out[0], fa);
else print('V', out[0].join("\t"));
out.shift();
}
//warn(tot_len + " alignment columns considered in calling");
warn(c1_len + " reference bases covered by exactly one contig");
warn(n_sub[0] + " substitutions; ts/tv = " + (n_sub[1]/n_sub[2]).toFixed(3));
warn(n_del[0] + " 1bp deletions");
warn(n_ins[0] + " 1bp insertions");
warn(n_del[1] + " 2bp deletions");
warn(n_ins[1] + " 2bp insertions");
warn(n_del[2] + " [3,"+gap_thres+") deletions");
warn(n_ins[2] + " [3,"+gap_thres+") insertions");
warn(n_del[3] + " ["+gap_thres+","+gap_thres_long+") deletions");
warn(n_ins[3] + " ["+gap_thres+","+gap_thres_long+") insertions");
warn(n_del[4] + " >=" + gap_thres_long + " deletions");
warn(n_ins[4] + " >=" + gap_thres_long + " insertions");
buf.destroy();
file.close();
if (fa != null) fasta_free(fa);
}
function paf_asmstat(args)
{
var c, min_query_len = 0, min_seg_len = 10000, max_diff = 0.01, bp_flank_len = 0, bp_gap_len = 0;
while ((c = getopt(args, "l:d:b:g:q:")) != null) {
if (c == 'l') min_seg_len = parseInt(getopt.arg);
else if (c == 'd') max_diff = parseFloat(getopt.arg);
else if (c == 'b') bp_flank_len = parseInt(getopt.arg);
else if (c == 'g') bp_gap_len = parseInt(getopt.arg);
else if (c == 'q') min_query_len = parseInt(getopt.arg);
}
if (getopt.ind == args.length) {
print("Usage: paftools.js asmstat [options] <ref.fa.fai> <asm1.paf> [...]");
print("Options:");
print(" -q INT ignore query shorter than INT [0]");
print(" -l INT min alignment block length [" + min_seg_len + "]");
print(" -d FLOAT max gap-compressed sequence divergence [" + max_diff + "]");
exit(1);
}
var file, buf = new Bytes();
var ref_len = 0;
file = new File(args[getopt.ind]);
while (file.readline(buf) >= 0) {
var t = buf.toString().split("\t");
ref_len += parseInt(t[1]);
}
file.close();
function process_query(qblocks, qblock_len, bp, qi) {
qblocks.sort(function(a,b) { return a[0]-b[0]; });
var last_k = null, last_blen = null, st = -1, en = -1, qcov = 0;
for (var k = 0; k < qblocks.length; ++k) {
var blen = qblocks[k][1] - qblocks[k][0];
if (k > 0 && qblocks[k][0] < qblocks[k-1][1]) {
if (qblocks[k][1] < qblocks[k-1][1]) continue;
blen = qblocks[k][1] - qblocks[k-1][1];
}
qblock_len.push(blen);
if (qblocks[k][0] > en) {
qcov += en - st;
st = qblocks[k][0];
en = qblocks[k][1];
} else en = en > qblocks[k][1]? en : qblocks[k][1];
if (last_k != null) {
var gap = 1000000000;
if (qblocks[k][2] == qblocks[last_k][2] && qblocks[k][3] == qblocks[last_k][3]) { // same chr and strand
var g1 = qblocks[k][0] - qblocks[last_k][1];
var g2 = qblocks[k][2] == '+'? qblocks[k][4] - qblocks[last_k][5] : qblocks[last_k][4] - qblocks[k][5];
gap = g1 > g2? g1 - g2 : g2 - g1;
}
var min = blen < last_blen? blen : last_blen;
var flank = k == 0? min : blen;
bp.push([flank, gap]);
qi.bp.push([flank, gap]);
}
last_k = k, last_blen = blen;
}
qcov += en - st;
return qcov;
}
function N50(lens, tot, quantile) {
lens.sort(function(a,b) { return b - a; });
if (tot == null) {
tot = 0;
for (var k = 0; k < lens.length; ++k)
tot += lens[k];
}
var sum = 0;
for (var k = 0; k < lens.length; ++k) {
if (sum <= quantile * tot && sum + lens[k] > quantile * tot)
return lens[k];
sum += lens[k];
}
}
function AUN(lens, tot) {
lens.sort(function(a,b) { return b - a; });
if (tot == null) {
tot = 0;
for (var k = 0; k < lens.length; ++k)
tot += lens[k];
}
var x = 0, y = 0;
for (var k = 0; k < lens.length; ++k) {
var l = x + lens[k] <= tot? lens[k] : tot - x;
x += lens[k];
y += l * (l / tot);
if (x >= tot) break;
}
return y.toFixed(0);
}
function count_bp(bp, min_blen, min_gap) {
var n_bp = 0;
for (var k = 0; k < bp.length; ++k)
if (bp[k][0] >= min_blen && bp[k][1] >= min_gap)
++n_bp;
return n_bp;
}
function compute_diff(cigar, NM) {
var m, re = /(\d+)([MID])/g;
var n_M = 0, n_gapo = 0, n_gaps = 0;
while ((m = re.exec(cigar)) != null) {
var len = parseInt(m[1]);
if (m[2] == 'M') n_M += len;
else ++n_gapo, n_gaps += len;
}
if (NM < n_gaps) throw Error('NM is smaller the number of gaps');
return (NM - n_gaps + n_gapo) / (n_M + n_gapo);
}
var labels = ['Length', 'l_cov', 'Rcov', 'Rdup', 'Qcov', 'NG75', 'NG50', 'NGA50', 'AUNGA', '#breaks', 'bp(' + min_seg_len + ',0)', 'bp(' + min_seg_len + ',10k)'];
var rst = [];
for (var i = 0; i < labels.length; ++i)
rst[i] = [];
var n_asm = args.length - (getopt.ind + 1);
var header = ["Metric"];
for (var i = 0; i < n_asm; ++i) {
var n_breaks = 0, qcov = 0;
var fn = args[getopt.ind + 1 + i];
var label = fn.replace(/.paf(.gz)?$/, "");
header.push(label);
var ref_blocks = [], qblock_len = [], qblocks = [], bp = [];
var query = {}, qinfo = {};
var last_qname = null;
file = new File(fn);
while (file.readline(buf) >= 0) {
var m, line = buf.toString();
var t = line.split("\t");
t[1] = parseInt(t[1]);
if (t[1] < min_query_len) continue;
if (t.length < 2) continue;
query[t[0]] = t[1];
if (qinfo[t[0]] == null) qinfo[t[0]] = {};
qinfo[t[0]].len = t[1];
qinfo[t[0]].bp = [];
if (t.length < 9 || t[5] == "*") continue;
if (!/\ttp:A:[PI]/.test(line)) continue;
var cigar = (m = /\tcg:Z:(\S+)/.exec(line)) != null? m[1] : null;
var NM = (m = /\tNM:i:(\d+)/.exec(line)) != null? parseInt(m[1]) : null;
var diff = cigar != null && NM != null? compute_diff(cigar, NM) : 0;
t[2] = parseInt(t[2]);
t[3] = parseInt(t[3]);
t[7] = parseInt(t[7]);
t[8] = parseInt(t[8]);
if (t[0] == last_qname) ++n_breaks;
if (diff > max_diff) continue;
if (t[3] - t[2] < min_seg_len) continue;
if (t[0] != last_qname) {
if (last_qname != null)
qcov += process_query(qblocks, qblock_len, bp, qinfo[last_qname]);
qblocks = [];
last_qname = t[0];
}
ref_blocks.push([t[5], t[7], t[8]]);
qblocks.push([t[2], t[3], t[4], t[5], t[7], t[8]]);
}
if (last_qname != null)
qcov += process_query(qblocks, qblock_len, bp, qinfo[last_qname]);
file.close();
// compute NG50
var asm_len = 0, asm_lens = []
for (var ctg in query) {
asm_len += query[ctg];
asm_lens.push(query[ctg]);
}
rst[0][i] = asm_len;
rst[5][i] = N50(asm_lens, ref_len, 0.75);
rst[6][i] = N50(asm_lens, ref_len, 0.5);
// compute coverage
var l_cov = 0;
ref_blocks.sort(function(a, b) { return a[0] > b[0]? 1 : a[0] < b[0]? -1 : a[1] - b[1]; });
var last_ref = null, st = -1, en = -1;
for (var j = 0; j < ref_blocks.length; ++j) {
if (ref_blocks[j][0] != last_ref || ref_blocks[j][1] > en) {
l_cov += en - st;
last_ref = ref_blocks[j][0];
st = ref_blocks[j][1];
en = ref_blocks[j][2];
} else en = en > ref_blocks[j][2]? en : ref_blocks[j][2];
}
l_cov += en - st;
rst[1][i] = l_cov;
rst[2][i] = (100.0 * (l_cov / ref_len)).toFixed(2) + '%';
rst[4][i] = (100.0 * (qcov / asm_len)).toFixed(2) + '%';
// compute cov1 and cov2+ lengths; see paf_call() for details
var c1_ctg = null, c1_start = 0, c1_end = 0, c1_len = 0;
for (var j = 0; j < ref_blocks.length; ++j) {
if (ref_blocks[j][0] != c1_ctg || ref_blocks[j][1] >= c1_end) {
if (c1_end > c1_start)
c1_len += c1_end - c1_start;
c1_ctg = ref_blocks[j][0], c1_start = ref_blocks[j][1], c1_end = ref_blocks[j][2];
} else if (ref_blocks[j][2] > c1_end) { // overlap
if (ref_blocks[j][1] > c1_start)
c1_len += ref_blocks[j][1] - c1_start;
c1_start = c1_end, c1_end = ref_blocks[j][2];
} else if (ref_blocks[j][2] > c1_start) { // contained
if (ref_blocks[j][1] > c1_start)
c1_len += ref_blocks[j][1] - c1_start;
c1_start = ref_blocks[j][2];
}
//print(ref_blocks[j][0], ref_blocks[j][1], ref_blocks[j][2], c1_start, c1_end, c1_len);
}
if (c1_end > c1_start)
c1_len += c1_end - c1_start;
rst[3][i] = (100 * (l_cov - c1_len) / l_cov).toFixed(2) + '%';
// compute NGA50
rst[7][i] = N50(qblock_len, ref_len, 0.5);
// compute AUNGA
rst[8][i] = AUN(qblock_len, ref_len);
// compute break points
rst[9][i] = n_breaks;
rst[10][i] = count_bp(bp, 500, 0);
rst[11][i] = count_bp(bp, 500, 10000);
// nb-plot; NOT USED
/*
var qa = [];
for (var qn in qinfo)
qa.push([qinfo[qn].len, qinfo[qn].bp]);
qa = qa.sort(function(a, b) { return b[0] - a[0] });
var sum = 0, n_bp = 0, next_quantile = 0.1;
for (var j = 0; j < qa.length; ++j) {
sum += qa[j][0];
for (var k = 0; k < qa[j][1].length; ++k)
if (qa[j][1][k][0] >= bp_flank_len && qa[j][1][k][1] >= bp_gap_len)
++n_bp;
if (sum >= ref_len * next_quantile) {
print(label, Math.floor(next_quantile * 100 + .5), qa[j][0], (sum / n_bp).toFixed(0), n_bp);
next_quantile += 0.1;
if (next_quantile >= 1.0) break;
}
}
*/
}
buf.destroy();
if (bp_flank_len <= 0) {
print(header.join("\t"));
for (var i = 0; i < labels.length; ++i)
print(labels[i], rst[i].join("\t"));
}
}
function paf_asmgene(args)
{
var c, opt = { min_cov:0.99, min_iden:0.99 }, print_err = false, auto_only = false;
while ((c = getopt(args, "i:c:ea")) != null)
if (c == 'i') opt.min_iden = parseFloat(getopt.arg);
else if (c == 'c') opt.min_cov = parseFloat(getopt.arg);
else if (c == 'e') print_err = true;
else if (c == 'a') auto_only = true;
var n_fn = args.length - getopt.ind;
if (n_fn < 2) {
print("Usage: paftools.js asmgene [options] <ref-splice.paf> <asm-splice.paf> [...]");
print("Options:");
print(" -i FLOAT min identity [" + opt.min_iden + "]");
print(" -c FLOAT min coverage [" + opt.min_cov + "]");
print(" -a only evaluate genes mapped to the autosomes");
print(" -e print fragmented/missing genes");
exit(1);
}
function process_query(opt, a) {
var b = [], cnt = [0, 0, 0];
for (var j = 0; j < a.length; ++j) {
if (a[j][4] < a[j][5] * opt.min_iden)
continue;
b.push(a[j].slice(0));
}
if (b.length == 0) return cnt;
// count full
var n_full = 0;
for (var j = 0; j < b.length; ++j)
if (b[j][3] - b[j][2] >= b[j][1] * opt.min_cov)
++n_full;
cnt[0] = n_full;
// compute coverage
b = b.sort(function(x, y) { return x[2] - y[2] });
var l_cov = 0, st = b[0][2], en = b[0][3];
for (var j = 1; j < b.length; ++j) {
if (b[j][2] <= en)
en = b[j][3] > en? b[j][3] : en;
else l_cov += en - st;
}
l_cov += en - st;
cnt[1] = l_cov / b[0][1];
cnt[2] = b.length;
return cnt;
}
var buf = new Bytes();
var gene = {}, header = [], refpos = {};
for (var i = getopt.ind; i < args.length; ++i) {
var fn = args[i];
var label = fn.replace(/.paf(.gz)?$/, "");
header.push(label);
var file = new File(fn), a = [];
while (file.readline(buf) >= 0) {
var t = buf.toString().split("\t");
var ql = parseInt(t[1]), qs = parseInt(t[2]), qe = parseInt(t[3]), mlen = parseInt(t[9]), blen = parseInt(t[10]), mapq = parseInt(t[11]);
if (i == getopt.ind) refpos[t[0]] = [t[0], t[1], t[5], t[7], t[8]];
if (gene[t[0]] == null) gene[t[0]] = [];
if (a.length && t[0] != a[0][0]) {
gene[a[0][0]][i - getopt.ind] = process_query(opt, a);
a = [];
}
a.push([t[0], ql, qs, qe, mlen, blen]);
}
if (a.length)
gene[t[0]][i - getopt.ind] = process_query(opt, a);
file.close();
}
// select the longest genes (not optimal, but should be good enough)
var gene_list = [], gene_nr = {};
for (var g in refpos)
gene_list.push(refpos[g]);
gene_list = gene_list.sort(function(a, b) { return a[2] < b[2]? -1 : a[2] > b[2]? 1 : a[3] - b[3] });
var last = 0;
for (var j = 1; j < gene_list.length; ++j) {
if (gene_list[j][2] != gene_list[last][2] || gene_list[j][3] >= gene_list[last][4]) {
gene_nr[gene_list[last][0]] = 1;
last = j;
} else if (gene_list[j][1] > gene_list[last][1]) {
last = j;
}
}
gene_nr[gene_list[last][0]] = 1;
// count and print
var col1 = ["full_sgl", "full_dup", "frag", "part50+", "part10+", "part10-", "dup_cnt", "dup_sum"];
var rst = [];
for (var k = 0; k < col1.length; ++k) {
rst[k] = [];
for (var i = 0; i < n_fn; ++i)
rst[k][i] = 0;
}
for (var g in gene) { // count single-copy genes
if (gene[g][0] == null || gene[g][0][0] != 1) continue;
if (gene_nr[g] == null) continue;
if (auto_only && /^(chr)?[XY]$/.test(refpos[g][2])) continue;
for (var i = 0; i < n_fn; ++i) {
if (gene[g][i] == null) {
rst[5][i]++;
if (print_err) print('M', header[i], refpos[g].join("\t"));
} else if (gene[g][i][0] == 1) {
rst[0][i]++;
} else if (gene[g][i][0] > 1) {
rst[1][i]++;
if (print_err) print('D', header[i], refpos[g].join("\t"));
} else if (gene[g][i][1] >= opt.min_cov) {
rst[2][i]++;
if (print_err) print('F', header[i], refpos[g].join("\t"));
} else if (gene[g][i][1] >= 0.5) {
rst[3][i]++;
if (print_err) print('5', header[i], refpos[g].join("\t"));
} else if (gene[g][i][1] >= 0.1) {
rst[4][i]++;
if (print_err) print('1', header[i], refpos[g].join("\t"));
} else {
rst[5][i]++;
if (print_err) print('0', header[i], refpos[g].join("\t")); // TODO: reduce code duplicates...
}
}
}
for (var g in gene) { // count multi-copy genes
if (gene[g][0] == null || gene[g][0][0] <= 1) continue;
if (gene_nr[g] == null) continue;
if (auto_only && /^(chr)?[XY]$/.test(refpos[g][2])) continue;
for (var i = 0; i < n_fn; ++i) {
if (gene[g][i] != null) rst[7][i] += gene[g][i][0];
if (gene[g][i] != null && gene[g][i][0] > 1) {
rst[6][i]++;
} else if (print_err) {
print('d', header[i], gene[g][0][0], refpos[g].join("\t"));
}
}
}
print('H', 'Metric', header.join("\t"));
for (var k = 0; k < rst.length; ++k) {
print('X', col1[k], rst[k].join("\t"));
}
buf.destroy();
}
function paf_stat(args)
{
var c, gap_out_len = null;
while ((c = getopt(args, "l:")) != null)
if (c == 'l') gap_out_len = parseInt(getopt.arg);
if (getopt.ind == args.length) {
print("Usage: paftools.js stat [-l gapOutLen] <in.sam>|<in.paf>");
exit(1);
}
var buf = new Bytes();
var file = args[getopt.ind] == '-'? new File() : new File(args[getopt.ind]);
var re = /(\d+)([MIDSHNX=])/g;
var lineno = 0, n_pri = 0, n_2nd = 0, n_seq = 0, n_cigar_64k = 0, l_tot = 0, l_cov = 0;
var n_gap = [[0, 0, 0, 0, 0, 0], [0, 0, 0, 0, 0, 0]];
function cov_len(regs)
{
regs.sort(function(a,b) {return a[0]-b[0]});
var st = regs[0][0], en = regs[0][1], l = 0;
for (var i = 1; i < regs.length; ++i) {
if (regs[i][0] < en)
en = en > regs[i][1]? en : regs[i][1];
else l += en - st, st = regs[i][0], en = regs[i][1];
}
l += en - st;
return l;
}
var last = null, last_qlen = null, regs = [];
while (file.readline(buf) >= 0) {
var line = buf.toString();
++lineno;
if (line.charAt(0) != '@') {
var t = line.split("\t", 12);
var m, rs, cigar = null, is_pri = false, is_sam = false, is_rev = false, tname = null;