Update debian changelog
[pysam.git] / samtools / phase.c.pysam.c
1 #include "pysam.h"
2
3 #include <stdio.h>
4 #include <stdlib.h>
5 #include <unistd.h>
6 #include <stdint.h>
7 #include <math.h>
8 #include <zlib.h>
9 #include "bam.h"
10 #include "errmod.h"
11
12 #include "kseq.h"
13 KSTREAM_INIT(gzFile, gzread, 16384)
14
15 #define MAX_VARS 256
16 #define FLIP_PENALTY 2
17 #define FLIP_THRES 4
18 #define MASK_THRES 3
19
20 #define FLAG_FIX_CHIMERA 0x1
21 #define FLAG_LIST_EXCL   0x4
22 #define FLAG_DROP_AMBI   0x8
23
24 typedef struct {
25         // configurations, initialized in the main function
26         int flag, k, min_baseQ, min_varLOD, max_depth;
27         // other global variables
28         int vpos_shift;
29         bamFile fp;
30         char *pre;
31         bamFile out[3];
32         // alignment queue
33         int n, m;
34         bam1_t **b;
35 } phaseg_t;
36
37 typedef struct {
38         int8_t seq[MAX_VARS]; // TODO: change to dynamic memory allocation!
39         int vpos, beg, end;
40         uint32_t vlen:16, single:1, flip:1, phase:1, phased:1, ambig:1;
41         uint32_t in:16, out:16; // in-phase and out-phase
42 } frag_t, *frag_p;
43
44 #define rseq_lt(a,b) ((a)->vpos < (b)->vpos)
45
46 #include "khash.h"
47 KHASH_SET_INIT_INT64(set64)
48 KHASH_MAP_INIT_INT64(64, frag_t)
49
50 typedef khash_t(64) nseq_t;
51
52 #include "ksort.h"
53 KSORT_INIT(rseq, frag_p, rseq_lt)
54
55 static char nt16_nt4_table[] = { 4, 0, 1, 4, 2, 4, 4, 4, 3, 4, 4, 4, 4, 4, 4, 4 };
56
57 static inline uint64_t X31_hash_string(const char *s)
58 {
59         uint64_t h = *s;
60         if (h) for (++s ; *s; ++s) h = (h << 5) - h + *s;
61         return h;
62 }
63
64 static void count1(int l, const uint8_t *seq, int *cnt)
65 {
66         int i, j, n_ambi;
67         uint32_t z, x;
68         if (seq[l-1] == 0) return; // do nothing is the last base is ambiguous
69         for (i = n_ambi = 0; i < l; ++i) // collect ambiguous bases
70                 if (seq[i] == 0) ++n_ambi;
71         if (l - n_ambi <= 1) return; // only one SNP
72         for (x = 0; x < 1u<<n_ambi; ++x) { // count
73                 for (i = j = 0, z = 0; i < l; ++i) {
74                         int c;
75                         if (seq[i]) c = seq[i] - 1;
76                         else {
77                                 c = x>>j&1;
78                                 ++j;
79                         }
80                         z = z<<1 | c;
81                 }
82                 ++cnt[z];
83         }
84 }
85
86 static int **count_all(int l, int vpos, nseq_t *hash)
87 {
88         khint_t k;
89         int i, j, **cnt;
90         uint8_t *seq;
91         seq = calloc(l, 1);
92         cnt = calloc(vpos, sizeof(void*));
93         for (i = 0; i < vpos; ++i) cnt[i] = calloc(1<<l, sizeof(int));
94         for (k = 0; k < kh_end(hash); ++k) {
95                 if (kh_exist(hash, k)) {
96                         frag_t *f = &kh_val(hash, k);
97                         if (f->vpos >= vpos || f->single) continue; // out of region; or singleton
98                         if (f->vlen == 1) { // such reads should be flagged as deleted previously if everything is right
99                                 f->single = 1;
100                                 continue;
101                         }
102                         for (j = 1; j < f->vlen; ++j) {
103                                 for (i = 0; i < l; ++i)
104                                         seq[i] = j < l - 1 - i? 0 : f->seq[j - (l - 1 - i)];
105                                 count1(l, seq, cnt[f->vpos + j]);
106                         }
107                 }
108         }
109         free(seq);
110         return cnt;
111 }
112
113 // phasing
114 static int8_t *dynaprog(int l, int vpos, int **w)
115 {
116         int *f[2], *curr, *prev, max, i;
117         int8_t **b, *h = 0;
118         uint32_t x, z = 1u<<(l-1), mask = (1u<<l) - 1;
119         f[0] = calloc(z, sizeof(int));
120         f[1] = calloc(z, sizeof(int));
121         b = calloc(vpos, sizeof(void*));
122         prev = f[0]; curr = f[1];
123         // fill the backtrack matrix
124         for (i = 0; i < vpos; ++i) {
125                 int *wi = w[i], *tmp;
126                 int8_t *bi;
127                 bi = b[i] = calloc(z, 1);
128                 /* In the following, x is the current state, which is the
129                  * lexicographically smaller local haplotype. xc is the complement of
130                  * x, or the larger local haplotype; y0 and y1 are the two predecessors
131                  * of x. */
132                 for (x = 0; x < z; ++x) { // x0 is the smaller 
133                         uint32_t y0, y1, xc;
134                         int c0, c1;
135                         xc = ~x&mask; y0 = x>>1; y1 = xc>>1;
136                         c0 = prev[y0] + wi[x] + wi[xc];
137                         c1 = prev[y1] + wi[x] + wi[xc];
138                         if (c0 > c1) bi[x] = 0, curr[x] = c0;
139                         else bi[x] = 1, curr[x] = c1;
140                 }
141                 tmp = prev; prev = curr; curr = tmp; // swap
142         }
143         { // backtrack
144                 uint32_t max_x = 0;
145                 int which = 0;
146                 h = calloc(vpos, 1);
147                 for (x = 0, max = 0, max_x = 0; x < z; ++x)
148                         if (prev[x] > max) max = prev[x], max_x = x;
149                 for (i = vpos - 1, x = max_x; i >= 0; --i) {
150                         h[i] = which? (~x&1) : (x&1);
151                         which = b[i][x]? !which : which;
152                         x = b[i][x]? (~x&mask)>>1 : x>>1;
153                 }
154         }
155         // free
156         for (i = 0; i < vpos; ++i) free(b[i]);
157         free(f[0]); free(f[1]); free(b);
158         return h;
159 }
160
161 // phase each fragment
162 static uint64_t *fragphase(int vpos, const int8_t *path, nseq_t *hash, int flip)
163 {
164         khint_t k;
165         uint64_t *pcnt;
166         uint32_t *left, *rght, max;
167         left = rght = 0; max = 0;
168         pcnt = calloc(vpos, 8);
169         for (k = 0; k < kh_end(hash); ++k) {
170                 if (kh_exist(hash, k)) {
171                         int i, c[2];
172                         frag_t *f = &kh_val(hash, k);
173                         if (f->vpos >= vpos) continue;
174                         // get the phase
175                         c[0] = c[1] = 0;
176                         for (i = 0; i < f->vlen; ++i) {
177                                 if (f->seq[i] == 0) continue;
178                                 ++c[f->seq[i] == path[f->vpos + i] + 1? 0 : 1];
179                         }
180                         f->phase = c[0] > c[1]? 0 : 1;
181                         f->in = c[f->phase]; f->out = c[1 - f->phase];
182                         f->phased = f->in == f->out? 0 : 1;
183                         f->ambig = (f->in && f->out && f->out < 3 && f->in <= f->out + 1)? 1 : 0;
184                         // fix chimera
185                         f->flip = 0;
186                         if (flip && c[0] >= 3 && c[1] >= 3) {
187                                 int sum[2], m, mi, md;
188                                 if (f->vlen > max) { // enlarge the array
189                                         max = f->vlen;
190                                         kroundup32(max);
191                                         left = realloc(left, max * 4);
192                                         rght = realloc(rght, max * 4);
193                                 }
194                                 for (i = 0, sum[0] = sum[1] = 0; i < f->vlen; ++i) { // get left counts
195                                         if (f->seq[i]) {
196                                                 int c = f->phase? 2 - f->seq[i] : f->seq[i] - 1;
197                                                 ++sum[c == path[f->vpos + i]? 0 : 1];
198                                         }
199                                         left[i] = sum[1]<<16 | sum[0];
200                                 }
201                                 for (i = f->vlen - 1, sum[0] = sum[1] = 0; i >= 0; --i) { // get right counts
202                                         if (f->seq[i]) {
203                                                 int c = f->phase? 2 - f->seq[i] : f->seq[i] - 1;
204                                                 ++sum[c == path[f->vpos + i]? 0 : 1];
205                                         }
206                                         rght[i] = sum[1]<<16 | sum[0];
207                                 }
208                                 // find the best flip point
209                                 for (i = m = 0, mi = -1, md = -1; i < f->vlen - 1; ++i) {
210                                         int a[2];
211                                         a[0] = (left[i]&0xffff) + (rght[i+1]>>16&0xffff) - (rght[i+1]&0xffff) * FLIP_PENALTY;
212                                         a[1] = (left[i]>>16&0xffff) + (rght[i+1]&0xffff) - (rght[i+1]>>16&0xffff) * FLIP_PENALTY;
213                                         if (a[0] > a[1]) {
214                                                 if (a[0] > m) m = a[0], md = 0, mi = i;
215                                         } else {
216                                                 if (a[1] > m) m = a[1], md = 1, mi = i;
217                                         }
218                                 }
219                                 if (m - c[0] >= FLIP_THRES && m - c[1] >= FLIP_THRES) { // then flip
220                                         f->flip = 1;
221                                         if (md == 0) { // flip the tail
222                                                 for (i = mi + 1; i < f->vlen; ++i)
223                                                         if (f->seq[i] == 1) f->seq[i] = 2;
224                                                         else if (f->seq[i] == 2) f->seq[i] = 1;
225                                         } else { // flip the head
226                                                 for (i = 0; i <= mi; ++i)
227                                                         if (f->seq[i] == 1) f->seq[i] = 2;
228                                                         else if (f->seq[i] == 2) f->seq[i] = 1;
229                                         }
230                                 }
231                         }
232                         // update pcnt[]
233                         if (!f->single) {
234                                 for (i = 0; i < f->vlen; ++i) {
235                                         int c;
236                                         if (f->seq[i] == 0) continue;
237                                         c = f->phase? 2 - f->seq[i] : f->seq[i] - 1;
238                                         if (c == path[f->vpos + i]) {
239                                                 if (f->phase == 0) ++pcnt[f->vpos + i];
240                                                 else pcnt[f->vpos + i] += 1ull<<32;
241                                         } else {
242                                                 if (f->phase == 0) pcnt[f->vpos + i] += 1<<16;
243                                                 else pcnt[f->vpos + i] += 1ull<<48;
244                                         }
245                                 }
246                         }
247                 }
248         }
249         free(left); free(rght);
250         return pcnt;
251 }
252
253 static uint64_t *genmask(int vpos, const uint64_t *pcnt, int *_n)
254 {
255         int i, max = 0, max_i = -1, m = 0, n = 0, beg = 0, score = 0;
256         uint64_t *list = 0;
257         for (i = 0; i < vpos; ++i) {
258                 uint64_t x = pcnt[i];
259                 int c[4], pre = score, s;
260                 c[0] = x&0xffff; c[1] = x>>16&0xffff; c[2] = x>>32&0xffff; c[3] = x>>48&0xffff;
261                 s = (c[1] + c[3] == 0)? -(c[0] + c[2]) : (c[1] + c[3] - 1);
262                 if (c[3] > c[2]) s += c[3] - c[2];
263                 if (c[1] > c[0]) s += c[1] - c[0];
264                 score += s;
265                 if (score < 0) score = 0;
266                 if (pre == 0 && score > 0) beg = i; // change from zero to non-zero
267                 if ((i == vpos - 1 || score == 0) && max >= MASK_THRES) {
268                         if (n == m) {
269                                 m = m? m<<1 : 4;
270                                 list = realloc(list, m * 8);
271                         }
272                         list[n++] = (uint64_t)beg<<32 | max_i;
273                         i = max_i; // reset i to max_i
274                         score = 0;
275                 } else if (score > max) max = score, max_i = i;
276                 if (score == 0) max = 0;
277         }
278         *_n = n;
279         return list;
280 }
281
282 // trim heading and tailing ambiguous bases; mark deleted and remove sequence
283 static int clean_seqs(int vpos, nseq_t *hash)
284 {
285         khint_t k;
286         int ret = 0;
287         for (k = 0; k < kh_end(hash); ++k) {
288                 if (kh_exist(hash, k)) {
289                         frag_t *f = &kh_val(hash, k);
290                         int beg, end, i;
291                         if (f->vpos >= vpos) {
292                                 ret = 1;
293                                 continue;
294                         }
295                         for (i = 0; i < f->vlen; ++i)
296                                 if (f->seq[i] != 0) break;
297                         beg = i;
298                         for (i = f->vlen - 1; i >= 0; --i)
299                                 if (f->seq[i] != 0) break;
300                         end = i + 1;
301                         if (end - beg <= 0) kh_del(64, hash, k);
302                         else {
303                                 if (beg != 0) memmove(f->seq, f->seq + beg, end - beg);
304                                 f->vpos += beg; f->vlen = end - beg;
305                                 f->single = f->vlen == 1? 1 : 0;
306                         }
307                 }
308         }
309         return ret;
310 }
311
312 static void dump_aln(phaseg_t *g, int min_pos, const nseq_t *hash)
313 {
314         int i, is_flip, drop_ambi;
315         drop_ambi = g->flag & FLAG_DROP_AMBI;
316         is_flip = (drand48() < 0.5);
317         for (i = 0; i < g->n; ++i) {
318                 int end, which;
319                 uint64_t key;
320                 khint_t k;
321                 bam1_t *b = g->b[i];
322                 key = X31_hash_string(bam1_qname(b));
323                 end = bam_calend(&b->core, bam1_cigar(b));
324                 if (end > min_pos) break;
325                 k = kh_get(64, hash, key);
326                 if (k == kh_end(hash)) which = 3;
327                 else {
328                         frag_t *f = &kh_val(hash, k);
329                         if (f->ambig) which = drop_ambi? 2 : 3;
330                         else if (f->phased && f->flip) which = 2;
331                         else if (f->phased == 0) which = 3;
332                         else { // phased and not flipped
333                                 char c = 'Y';
334                                 which = f->phase;
335                                 bam_aux_append(b, "ZP", 'A', 1, (uint8_t*)&c);
336                         }
337                         if (which < 2 && is_flip) which = 1 - which; // increase the randomness
338                 }
339                 if (which == 3) which = (drand48() < 0.5);
340                 bam_write1(g->out[which], b);
341                 bam_destroy1(b);
342                 g->b[i] = 0;
343         }
344         memmove(g->b, g->b + i, (g->n - i) * sizeof(void*));
345         g->n -= i;
346 }
347
348 static int phase(phaseg_t *g, const char *chr, int vpos, uint64_t *cns, nseq_t *hash)
349 {
350         int i, j, n_seqs = kh_size(hash), n_masked = 0, min_pos;
351         khint_t k;
352         frag_t **seqs;
353         int8_t *path, *sitemask;
354         uint64_t *pcnt, *regmask;
355
356         if (vpos == 0) return 0;
357         i = clean_seqs(vpos, hash); // i is true if hash has an element with its vpos >= vpos
358         min_pos = i? cns[vpos]>>32 : 0x7fffffff;
359         if (vpos == 1) {
360                 printf("PS\t%s\t%d\t%d\n", chr, (int)(cns[0]>>32) + 1, (int)(cns[0]>>32) + 1);
361                 printf("M0\t%s\t%d\t%d\t%c\t%c\t%d\t0\t0\t0\t0\n//\n", chr, (int)(cns[0]>>32) + 1, (int)(cns[0]>>32) + 1,
362                         "ACGTX"[cns[0]&3], "ACGTX"[cns[0]>>16&3], g->vpos_shift + 1);
363                 for (k = 0; k < kh_end(hash); ++k) {
364                         if (kh_exist(hash, k)) {
365                                 frag_t *f = &kh_val(hash, k);
366                                 if (f->vpos) continue;
367                                 f->flip = 0;
368                                 if (f->seq[0] == 0) f->phased = 0;
369                                 else f->phased = 1, f->phase = f->seq[0] - 1;
370                         }
371                 }
372                 dump_aln(g, min_pos, hash);
373                 ++g->vpos_shift;
374                 return 1;
375         }
376         { // phase
377                 int **cnt;
378                 uint64_t *mask;
379                 printf("PS\t%s\t%d\t%d\n", chr, (int)(cns[0]>>32) + 1, (int)(cns[vpos-1]>>32) + 1);
380                 sitemask = calloc(vpos, 1);
381                 cnt = count_all(g->k, vpos, hash);
382                 path = dynaprog(g->k, vpos, cnt);
383                 for (i = 0; i < vpos; ++i) free(cnt[i]);
384                 free(cnt);
385                 pcnt = fragphase(vpos, path, hash, 0); // do not fix chimeras when masking
386                 mask = genmask(vpos, pcnt, &n_masked);
387                 regmask = calloc(n_masked, 8);
388                 for (i = 0; i < n_masked; ++i) {
389                         regmask[i] = cns[mask[i]>>32]>>32<<32 | cns[(uint32_t)mask[i]]>>32;
390                         for (j = mask[i]>>32; j <= (int32_t)mask[i]; ++j)
391                                 sitemask[j] = 1;
392                 }
393                 free(mask);
394                 if (g->flag & FLAG_FIX_CHIMERA) {
395                         free(pcnt);
396                         pcnt = fragphase(vpos, path, hash, 1);
397                 }
398         }
399         for (i = 0; i < n_masked; ++i)
400                 printf("FL\t%s\t%d\t%d\n", chr, (int)(regmask[i]>>32) + 1, (int)regmask[i] + 1);
401         for (i = 0; i < vpos; ++i) {
402                 uint64_t x = pcnt[i];
403                 int8_t c[2];
404                 c[0] = (cns[i]&0xffff)>>2 == 0? 4 : (cns[i]&3);
405                 c[1] = (cns[i]>>16&0xffff)>>2 == 0? 4 : (cns[i]>>16&3);
406                 printf("M%d\t%s\t%d\t%d\t%c\t%c\t%d\t%d\t%d\t%d\t%d\n", sitemask[i]+1, chr, (int)(cns[0]>>32) + 1, (int)(cns[i]>>32) + 1, "ACGTX"[c[path[i]]], "ACGTX"[c[1-path[i]]],
407                         i + g->vpos_shift + 1, (int)(x&0xffff), (int)(x>>16&0xffff), (int)(x>>32&0xffff), (int)(x>>48&0xffff));
408         }
409         free(path); free(pcnt); free(regmask); free(sitemask);
410         seqs = calloc(n_seqs, sizeof(void*));
411         for (k = 0, i = 0; k < kh_end(hash); ++k) 
412                 if (kh_exist(hash, k) && kh_val(hash, k).vpos < vpos && !kh_val(hash, k).single)
413                         seqs[i++] = &kh_val(hash, k);
414         n_seqs = i;
415         ks_introsort_rseq(n_seqs, seqs);
416         for (i = 0; i < n_seqs; ++i) {
417                 frag_t *f = seqs[i];
418                 printf("EV\t0\t%s\t%d\t40\t%dM\t*\t0\t0\t", chr, f->vpos + 1 + g->vpos_shift, f->vlen);
419                 for (j = 0; j < f->vlen; ++j) {
420                         uint32_t c = cns[f->vpos + j];
421                         if (f->seq[j] == 0) putchar('N');
422                         else putchar("ACGT"[f->seq[j] == 1? (c&3) : (c>>16&3)]);
423                 }
424                 printf("\t*\tYP:i:%d\tYF:i:%d\tYI:i:%d\tYO:i:%d\tYS:i:%d\n", f->phase, f->flip, f->in, f->out, f->beg+1);
425         }
426         free(seqs);
427         printf("//\n");
428         fflush(stdout);
429         g->vpos_shift += vpos;
430         dump_aln(g, min_pos, hash);
431         return vpos;
432 }
433
434 static void update_vpos(int vpos, nseq_t *hash)
435 {
436         khint_t k;
437         for (k = 0; k < kh_end(hash); ++k) {
438                 if (kh_exist(hash, k)) {
439                         frag_t *f = &kh_val(hash, k);
440                         if (f->vpos < vpos) kh_del(64, hash, k); // TODO: if frag_t::seq is allocated dynamically, free it
441                         else f->vpos -= vpos;
442                 }
443         }
444 }
445
446 static nseq_t *shrink_hash(nseq_t *hash) // TODO: to implement
447 {
448         return hash;
449 }
450
451 static int readaln(void *data, bam1_t *b)
452 {
453         phaseg_t *g = (phaseg_t*)data;
454         int ret;
455         ret = bam_read1(g->fp, b);
456         if (ret < 0) return ret;
457         if (!(b->core.flag & (BAM_FUNMAP|BAM_FSECONDARY|BAM_FQCFAIL|BAM_FDUP)) && g->pre) {
458                 if (g->n == g->m) {
459                         g->m = g->m? g->m<<1 : 16;
460                         g->b = realloc(g->b, g->m * sizeof(void*));
461                 }
462                 g->b[g->n++] = bam_dup1(b);
463         }
464         return ret;
465 }
466
467 static khash_t(set64) *loadpos(const char *fn, bam_header_t *h)
468 {
469         gzFile fp;
470         kstream_t *ks;
471         int ret, dret;
472         kstring_t *str;
473         khash_t(set64) *hash;
474
475         hash = kh_init(set64);
476         str = calloc(1, sizeof(kstring_t));
477         fp = strcmp(fn, "-")? gzopen(fn, "r") : gzdopen(fileno(stdin), "r");
478         ks = ks_init(fp);
479         while (ks_getuntil(ks, 0, str, &dret) >= 0) {
480                 int tid = bam_get_tid(h, str->s);
481                 if (tid >= 0 && dret != '\n') {
482                         if (ks_getuntil(ks, 0, str, &dret) >= 0) {
483                                 uint64_t x = (uint64_t)tid<<32 | (atoi(str->s) - 1);
484                                 kh_put(set64, hash, x, &ret);
485                         } else break;
486                 }
487                 if (dret != '\n') while ((dret = ks_getc(ks)) > 0 && dret != '\n');
488                 if (dret < 0) break;
489         }
490         ks_destroy(ks);
491         gzclose(fp);
492         free(str->s); free(str);
493         return hash;
494 }
495
496 static int gl2cns(float q[16])
497 {
498         int i, j, min_ij;
499         float min, min2;
500         min = min2 = 1e30; min_ij = -1;
501         for (i = 0; i < 4; ++i) {
502                 for (j = i; j < 4; ++j) {
503                         if (q[i<<2|j] < min) min_ij = i<<2|j, min2 = min, min = q[i<<2|j];
504                         else if (q[i<<2|j] < min2) min2 = q[i<<2|j];
505                 }
506         }
507         return (min_ij>>2&3) == (min_ij&3)? 0 : 1<<18 | (min_ij>>2&3)<<16 | (min_ij&3) | (int)(min2 - min + .499) << 2;
508 }
509
510 int main_phase(int argc, char *argv[])
511 {
512         extern void bam_init_header_hash(bam_header_t *header);
513         int c, tid, pos, vpos = 0, n, lasttid = -1, max_vpos = 0;
514         const bam_pileup1_t *plp;
515         bam_plp_t iter;
516         bam_header_t *h;
517         nseq_t *seqs;
518         uint64_t *cns = 0;
519         phaseg_t g;
520         char *fn_list = 0;
521         khash_t(set64) *set = 0;
522         errmod_t *em;
523         uint16_t *bases;
524
525         memset(&g, 0, sizeof(phaseg_t));
526         g.flag = FLAG_FIX_CHIMERA;
527         g.min_varLOD = 37; g.k = 13; g.min_baseQ = 13; g.max_depth = 256;
528         while ((c = getopt(argc, argv, "Q:eFq:k:b:l:D:A:")) >= 0) {
529                 switch (c) {
530                         case 'D': g.max_depth = atoi(optarg); break;
531                         case 'q': g.min_varLOD = atoi(optarg); break;
532                         case 'Q': g.min_baseQ = atoi(optarg); break;
533                         case 'k': g.k = atoi(optarg); break;
534                         case 'F': g.flag &= ~FLAG_FIX_CHIMERA; break;
535                         case 'e': g.flag |= FLAG_LIST_EXCL; break;
536                         case 'A': g.flag |= FLAG_DROP_AMBI; break;
537                         case 'b': g.pre = strdup(optarg); break;
538                         case 'l': fn_list = strdup(optarg); break;
539                 }
540         }
541         if (argc == optind) {
542                 fprintf(pysamerr, "\n");
543                 fprintf(pysamerr, "Usage:   samtools phase [options] <in.bam>\n\n");
544                 fprintf(pysamerr, "Options: -k INT    block length [%d]\n", g.k);
545                 fprintf(pysamerr, "         -b STR    prefix of BAMs to output [null]\n");
546                 fprintf(pysamerr, "         -q INT    min het phred-LOD [%d]\n", g.min_varLOD);
547                 fprintf(pysamerr, "         -Q INT    min base quality in het calling [%d]\n", g.min_baseQ);
548                 fprintf(pysamerr, "         -D INT    max read depth [%d]\n", g.max_depth);
549 //              fprintf(pysamerr, "         -l FILE   list of sites to phase [null]\n");
550                 fprintf(pysamerr, "         -F        do not attempt to fix chimeras\n");
551                 fprintf(pysamerr, "         -A        drop reads with ambiguous phase\n");
552 //              fprintf(pysamerr, "         -e        do not discover SNPs (effective with -l)\n");
553                 fprintf(pysamerr, "\n");
554                 return 1;
555         }
556         g.fp = strcmp(argv[optind], "-")? bam_open(argv[optind], "r") : bam_dopen(fileno(stdin), "r");
557         h = bam_header_read(g.fp);
558         if (fn_list) { // read the list of sites to phase
559                 bam_init_header_hash(h);
560                 set = loadpos(fn_list, h);
561                 free(fn_list);
562         } else g.flag &= ~FLAG_LIST_EXCL;
563         if (g.pre) { // open BAMs to write
564                 char *s = malloc(strlen(g.pre) + 20);
565                 strcpy(s, g.pre); strcat(s, ".0.bam"); g.out[0] = bam_open(s, "w");
566                 strcpy(s, g.pre); strcat(s, ".1.bam"); g.out[1] = bam_open(s, "w");
567                 strcpy(s, g.pre); strcat(s, ".chimera.bam"); g.out[2] = bam_open(s, "w");
568                 for (c = 0; c <= 2; ++c) bam_header_write(g.out[c], h);
569                 free(s);
570         }
571
572         iter = bam_plp_init(readaln, &g);
573         g.vpos_shift = 0;
574         seqs = kh_init(64);
575         em = errmod_init(1. - 0.83);
576         bases = calloc(g.max_depth, 2);
577         printf("CC\n");
578         printf("CC\tDescriptions:\nCC\n");
579         printf("CC\t  CC      comments\n");
580         printf("CC\t  PS      start of a phase set\n");
581         printf("CC\t  FL      filtered region\n");
582         printf("CC\t  M[012]  markers; 0 for singletons, 1 for phased and 2 for filtered\n");
583         printf("CC\t  EV      supporting reads; SAM format\n");
584         printf("CC\t  //      end of a phase set\nCC\n");
585         printf("CC\tFormats of PS, FL and M[012] lines (1-based coordinates):\nCC\n");
586         printf("CC\t  PS  chr  phaseSetStart  phaseSetEnd\n");
587         printf("CC\t  FL  chr  filterStart    filterEnd\n");
588         printf("CC\t  M?  chr  PS  pos  allele0  allele1  hetIndex  #supports0  #errors0  #supp1  #err1\n");
589         printf("CC\nCC\n");
590         fflush(stdout);
591         while ((plp = bam_plp_auto(iter, &tid, &pos, &n)) != 0) {
592                 int i, k, c, tmp, dophase = 1, in_set = 0;
593                 float q[16];
594                 if (tid < 0) break;
595                 if (tid != lasttid) { // change of chromosome
596                         g.vpos_shift = 0;
597                         if (lasttid >= 0) {
598                                 seqs = shrink_hash(seqs);
599                                 phase(&g, h->target_name[lasttid], vpos, cns, seqs);
600                                 update_vpos(0x7fffffff, seqs);
601                         }
602                         lasttid = tid;
603                         vpos = 0;
604                 }
605                 if (set && kh_get(set64, set, (uint64_t)tid<<32 | pos) != kh_end(set)) in_set = 1;
606                 if (n > g.max_depth) continue; // do not proceed if the depth is too high
607                 // fill the bases array and check if there is a variant
608                 for (i = k = 0; i < n; ++i) {
609                         const bam_pileup1_t *p = plp + i;
610                         uint8_t *seq;
611                         int q, baseQ, b;
612                         if (p->is_del || p->is_refskip) continue;
613                         baseQ = bam1_qual(p->b)[p->qpos];
614                         if (baseQ < g.min_baseQ) continue;
615                         seq = bam1_seq(p->b);
616                         b = bam_nt16_nt4_table[bam1_seqi(seq, p->qpos)];
617                         if (b > 3) continue;
618                         q = baseQ < p->b->core.qual? baseQ : p->b->core.qual;
619                         if (q < 4) q = 4;
620                         if (q > 63) q = 63;
621                         bases[k++] = q<<5 | (int)bam1_strand(p->b)<<4 | b;
622                 }
623                 if (k == 0) continue;
624                 errmod_cal(em, k, 4, bases, q); // compute genotype likelihood
625                 c = gl2cns(q); // get the consensus
626                 // tell if to proceed
627                 if (set && (g.flag&FLAG_LIST_EXCL) && !in_set) continue; // not in the list
628                 if (!in_set && (c&0xffff)>>2 < g.min_varLOD) continue; // not a variant
629                 // add the variant
630                 if (vpos == max_vpos) {
631                         max_vpos = max_vpos? max_vpos<<1 : 128;
632                         cns = realloc(cns, max_vpos * 8);
633                 }
634                 cns[vpos] = (uint64_t)pos<<32 | c;
635                 for (i = 0; i < n; ++i) {
636                         const bam_pileup1_t *p = plp + i;
637                         uint64_t key;
638                         khint_t k;
639                         uint8_t *seq = bam1_seq(p->b);
640                         frag_t *f;
641                         if (p->is_del || p->is_refskip) continue;
642                         if (p->b->core.qual == 0) continue;
643                         // get the base code
644                         c = nt16_nt4_table[(int)bam1_seqi(seq, p->qpos)];
645                         if (c == (cns[vpos]&3)) c = 1;
646                         else if (c == (cns[vpos]>>16&3)) c = 2;
647                         else c = 0;
648                         // write to seqs
649                         key = X31_hash_string(bam1_qname(p->b));
650                         k = kh_put(64, seqs, key, &tmp);
651                         f = &kh_val(seqs, k);
652                         if (tmp == 0) { // present in the hash table
653                                 if (vpos - f->vpos + 1 < MAX_VARS) {
654                                         f->vlen = vpos - f->vpos + 1;
655                                         f->seq[f->vlen-1] = c;
656                                         f->end = bam_calend(&p->b->core, bam1_cigar(p->b));
657                                 }
658                                 dophase = 0;
659                         } else { // absent
660                                 memset(f->seq, 0, MAX_VARS);
661                                 f->beg = p->b->core.pos;
662                                 f->end = bam_calend(&p->b->core, bam1_cigar(p->b));
663                                 f->vpos = vpos, f->vlen = 1, f->seq[0] = c, f->single = f->phased = f->flip = f->ambig = 0;
664                         }
665                 }
666                 if (dophase) {
667                         seqs = shrink_hash(seqs);
668                         phase(&g, h->target_name[tid], vpos, cns, seqs);
669                         update_vpos(vpos, seqs);
670                         cns[0] = cns[vpos];
671                         vpos = 0;
672                 }
673                 ++vpos;
674         }
675         if (tid >= 0) phase(&g, h->target_name[tid], vpos, cns, seqs);
676         bam_header_destroy(h);
677         bam_plp_destroy(iter);
678         bam_close(g.fp);
679         kh_destroy(64, seqs);
680         kh_destroy(set64, set);
681         free(cns);
682         errmod_destroy(em);
683         free(bases);
684         if (g.pre) {
685                 for (c = 0; c <= 2; ++c) bam_close(g.out[c]);
686                 free(g.pre); free(g.b);
687         }
688         return 0;
689 }