]> git.donarmstrong.com Git - samtools.git/blob - misc/bamcheck.c
Removed glibc dependency (getline)
[samtools.git] / misc / bamcheck.c
1 /* 
2     Author: petr.danecek@sanger
3     gcc -Wall -Winline -g -O2 -I ~/git/samtools bamcheck.c -o bamcheck -lm -lz -L ~/git/samtools -lbam
4
5     Assumptions, approximations and other issues:
6         - GC-depth graph does not split reads, the starting position determines which bin is incremented.
7             There are small overlaps between bins (max readlen-1). However, the bins are big (20k).
8         - coverage distribution ignores softclips and deletions
9         - some stats require sorted BAMs
10         - GC content graph can have an untidy, step-like pattern when BAM contains multiple read lengths.
11         - The whole reads are used with -t, no splicing is done, no indels or soft clips are 
12             considered, even small overlap is good enough to include the read in the stats.
13 */
14
15 #define BAMCHECK_VERSION "2012-04-04"
16
17 #define _ISOC99_SOURCE
18 #include <stdio.h>
19 #include <stdlib.h>
20 #include <stdarg.h>
21 #include <string.h>
22 #include <math.h>
23 #include <ctype.h>
24 #include <getopt.h>
25 #include <errno.h>
26 #include "faidx.h"
27 #include "khash.h"
28 #include "sam.h"
29 #include "razf.h"
30
31 #define BWA_MIN_RDLEN 35
32 #define IS_PAIRED(bam) ((bam)->core.flag&BAM_FPAIRED && !((bam)->core.flag&BAM_FUNMAP) && !((bam)->core.flag&BAM_FMUNMAP))
33 #define IS_UNMAPPED(bam) ((bam)->core.flag&BAM_FUNMAP)
34 #define IS_REVERSE(bam) ((bam)->core.flag&BAM_FREVERSE)
35 #define IS_MATE_REVERSE(bam) ((bam)->core.flag&BAM_FMREVERSE)
36 #define IS_READ1(bam) ((bam)->core.flag&BAM_FREAD1)
37 #define IS_READ2(bam) ((bam)->core.flag&BAM_FREAD2)
38 #define IS_DUP(bam) ((bam)->core.flag&BAM_FDUP)
39
40 typedef struct 
41 {
42     int32_t line_len, line_blen;
43     int64_t len;
44     uint64_t offset;
45
46 faidx1_t;
47 KHASH_MAP_INIT_STR(s, faidx1_t)
48 KHASH_MAP_INIT_STR(str, int)
49 struct __faidx_t {
50     RAZF *rz;
51     int n, m;
52     char **name;
53     khash_t(s) *hash;
54 };
55
56 typedef struct
57 {
58     float gc;
59     uint32_t depth;
60 }
61 gc_depth_t;
62
63 // For coverage distribution, a simple pileup
64 typedef struct 
65 {
66     int64_t pos;
67     int size, start;
68     int *buffer;
69
70 round_buffer_t;
71
72 typedef struct { uint32_t from, to; } pos_t;
73 typedef struct
74 {
75     int npos,mpos,cpos;
76     pos_t *pos;
77 }
78 regions_t;
79
80 typedef struct
81 {
82     // Parameters
83     int trim_qual;      // bwa trim quality
84     int rmdup;          // Exclude reads marked as duplicates from the stats
85
86     // Dimensions of the quality histogram holder (quals_1st,quals_2nd), GC content holder (gc_1st,gc_2nd),
87     //  insert size histogram holder
88     int nquals;         // The number of quality bins 
89     int nbases;         // The maximum sequence length the allocated array can hold
90     int nisize;         // The maximum insert size that the allocated array can hold
91     int ngc;            // The size of gc_1st and gc_2nd
92     int nindels;        // The maximum indel length for indel distribution
93
94     // Arrays for the histogram data
95     uint64_t *quals_1st, *quals_2nd;
96     uint64_t *gc_1st, *gc_2nd;
97     uint64_t *isize_inward, *isize_outward, *isize_other;
98     uint64_t *acgt_cycles;
99     uint64_t *read_lengths;
100     uint64_t *insertions, *deletions;
101     uint64_t *ins_cycles, *del_cycles;
102
103     // The extremes encountered
104     int max_len;            // Maximum read length
105     int max_qual;           // Maximum quality
106     float isize_main_bulk;  // There are always some unrealistically big insert sizes, report only the main part
107     int is_sorted;
108
109     // Summary numbers
110     uint64_t total_len;
111     uint64_t total_len_dup;
112     uint64_t nreads_1st;
113     uint64_t nreads_2nd;
114     uint64_t nreads_dup;
115     uint64_t nreads_unmapped;
116     uint64_t nreads_unpaired;
117     uint64_t nreads_paired;
118     uint64_t nreads_mq0;
119     uint64_t nbases_mapped;
120     uint64_t nbases_mapped_cigar;
121     uint64_t nbases_trimmed;  // bwa trimmed bases
122     uint64_t nmismatches;
123
124     // GC-depth related data
125     uint32_t ngcd, igcd;        // The maximum number of GC depth bins and index of the current bin
126     gc_depth_t *gcd;            // The GC-depth bins holder
127     int gcd_bin_size;           // The size of GC-depth bin
128     int32_t tid, gcd_pos;       // Position of the current bin
129     int32_t pos;                // Position of the last read
130
131     // Coverage distribution related data
132     int ncov;                       // The number of coverage bins
133     uint64_t *cov;                  // The coverage frequencies
134     int cov_min,cov_max,cov_step;   // Minimum, maximum coverage and size of the coverage bins
135     round_buffer_t cov_rbuf;        // Pileup round buffer
136
137     // Mismatches by read cycle 
138     uint8_t *rseq_buf;             // A buffer for reference sequence to check the mismatches against
139     int nref_seq;               // The size of the buffer
140     int32_t rseq_pos;           // The coordinate of the first base in the buffer
141     int32_t rseq_len;           // The used part of the buffer
142     uint64_t *mpc_buf;          // Mismatches per cycle
143
144     // Filters
145     int filter_readlen;
146
147     // Target regions
148     int nregions;
149     regions_t *regions;
150
151     // Auxiliary data
152     double sum_qual;            // For calculating average quality value 
153     samfile_t *sam;             // Unused
154     faidx_t *fai;               // Reference sequence for GC-depth graph
155     int argc;                   // Command line arguments to be printed on the output
156     char **argv;
157 }
158 stats_t;
159
160 void error(const char *format, ...);
161
162 // Coverage distribution methods
163 inline int coverage_idx(int min, int max, int n, int step, int depth)
164 {
165     if ( depth < min )
166         return 0;
167
168     if ( depth > max )
169         return n-1;
170
171     return 1 + (depth - min) / step;
172 }
173
174 inline int round_buffer_lidx2ridx(int offset, int size, int64_t refpos, int64_t pos)
175 {
176     return (offset + (pos-refpos) % size) % size;
177 }
178
179 void round_buffer_flush(stats_t *stats, int64_t pos)
180 {
181     int ibuf,idp;
182
183     if ( pos==stats->cov_rbuf.pos ) 
184         return;
185
186     int64_t new_pos = pos;
187     if ( pos==-1 || pos - stats->cov_rbuf.pos >= stats->cov_rbuf.size )
188     {
189         // Flush the whole buffer, but in sequential order, 
190         pos = stats->cov_rbuf.pos + stats->cov_rbuf.size - 1;
191     }
192
193     if ( pos < stats->cov_rbuf.pos ) 
194         error("Expected coordinates in ascending order, got %ld after %ld\n", pos,stats->cov_rbuf.pos);
195
196     int ifrom = stats->cov_rbuf.start;
197     int ito = round_buffer_lidx2ridx(stats->cov_rbuf.start,stats->cov_rbuf.size,stats->cov_rbuf.pos,pos-1);
198     if ( ifrom>ito )
199     {
200         for (ibuf=ifrom; ibuf<stats->cov_rbuf.size; ibuf++)
201         {
202             if ( !stats->cov_rbuf.buffer[ibuf] )
203                 continue;
204             idp = coverage_idx(stats->cov_min,stats->cov_max,stats->ncov,stats->cov_step,stats->cov_rbuf.buffer[ibuf]);
205             stats->cov[idp]++;
206             stats->cov_rbuf.buffer[ibuf] = 0;
207         }
208         ifrom = 0;
209     }
210     for (ibuf=ifrom; ibuf<=ito; ibuf++)
211     {
212         if ( !stats->cov_rbuf.buffer[ibuf] )
213             continue;
214         idp = coverage_idx(stats->cov_min,stats->cov_max,stats->ncov,stats->cov_step,stats->cov_rbuf.buffer[ibuf]);
215         stats->cov[idp]++;
216         stats->cov_rbuf.buffer[ibuf] = 0;
217     }
218     stats->cov_rbuf.start = (new_pos==-1) ? 0 : round_buffer_lidx2ridx(stats->cov_rbuf.start,stats->cov_rbuf.size,stats->cov_rbuf.pos,pos);
219     stats->cov_rbuf.pos   = new_pos;
220 }
221
222 void round_buffer_insert_read(round_buffer_t *rbuf, int64_t from, int64_t to)
223 {
224     if ( to-from >= rbuf->size )
225         error("The read length too big (%d), please increase the buffer length (currently %d)\n", to-from+1,rbuf->size);
226     if ( from < rbuf->pos )
227         error("The reads are not sorted (%ld comes after %ld).\n", from,rbuf->pos);
228
229     int ifrom,ito,ibuf;
230     ifrom = round_buffer_lidx2ridx(rbuf->start,rbuf->size,rbuf->pos,from);
231     ito   = round_buffer_lidx2ridx(rbuf->start,rbuf->size,rbuf->pos,to);
232     if ( ifrom>ito )
233     {
234         for (ibuf=ifrom; ibuf<rbuf->size; ibuf++)
235             rbuf->buffer[ibuf]++;
236         ifrom = 0;
237     }
238     for (ibuf=ifrom; ibuf<=ito; ibuf++)
239         rbuf->buffer[ibuf]++;
240 }
241
242 // Calculate the number of bases in the read trimmed by BWA
243 int bwa_trim_read(int trim_qual, uint8_t *quals, int len, int reverse) 
244 {
245     if ( len<BWA_MIN_RDLEN ) return 0;
246
247     // Although the name implies that the read cannot be trimmed to more than BWA_MIN_RDLEN,
248     //  the calculation can in fact trim it to (BWA_MIN_RDLEN-1). (bwa_trim_read in bwa/bwaseqio.c).
249     int max_trimmed = len - BWA_MIN_RDLEN + 1;
250     int l, sum=0, max_sum=0, max_l=0;
251
252     for (l=0; l<max_trimmed; l++)
253     {
254         sum += trim_qual - quals[ reverse ? l : len-1-l ];
255         if ( sum<0 ) break;
256         if ( sum>max_sum )
257         {
258             max_sum = sum;
259             // This is the correct way, but bwa clips from some reason one base less
260             // max_l   = l+1;
261             max_l   = l;
262         }
263     }
264     return max_l;
265 }
266
267
268 void count_indels(stats_t *stats,bam1_t *bam_line) 
269 {
270     int is_fwd = IS_REVERSE(bam_line) ? 0 : 1;
271     int icig;
272     int icycle = 0;
273     int read_len = bam_line->core.l_qseq;
274     for (icig=0; icig<bam_line->core.n_cigar; icig++) 
275     {
276         // Conversion from uint32_t to MIDNSHP
277         //  0123456
278         //  MIDNSHP
279         int cig  = bam1_cigar(bam_line)[icig] & BAM_CIGAR_MASK;
280         int ncig = bam1_cigar(bam_line)[icig] >> BAM_CIGAR_SHIFT;
281
282         if ( cig==1 )
283         {
284             int idx = is_fwd ? icycle : read_len-icycle-1;
285             if ( idx >= stats->nbases ) error("FIXME: %d vs %d\n", idx,stats->nbases);
286             stats->ins_cycles[idx]++;
287             icycle += ncig;
288             if ( ncig<=stats->nindels )
289                 stats->insertions[ncig-1]++;
290             continue;
291         }
292         if ( cig==2 )
293         {
294             int idx = is_fwd ? icycle : read_len-icycle-1;
295             if ( idx >= stats->nbases ) error("FIXME: %d vs %d\n", idx,stats->nbases);
296             stats->del_cycles[idx]++;
297             if ( ncig<=stats->nindels )
298                 stats->deletions[ncig-1]++;
299             continue;
300         }
301         icycle += ncig;
302     }
303 }
304
305 void count_mismatches_per_cycle(stats_t *stats,bam1_t *bam_line) 
306 {
307     int is_fwd = IS_REVERSE(bam_line) ? 0 : 1;
308     int icig,iread=0,icycle=0;
309     int iref = bam_line->core.pos - stats->rseq_pos;
310     int read_len   = bam_line->core.l_qseq;
311     uint8_t *read  = bam1_seq(bam_line);
312     uint8_t *quals = bam1_qual(bam_line);
313     uint64_t *mpc_buf = stats->mpc_buf;
314     for (icig=0; icig<bam_line->core.n_cigar; icig++) 
315     {
316         // Conversion from uint32_t to MIDNSHP
317         //  0123456
318         //  MIDNSHP
319         int cig  = bam1_cigar(bam_line)[icig] & BAM_CIGAR_MASK;
320         int ncig = bam1_cigar(bam_line)[icig] >> BAM_CIGAR_SHIFT;
321         if ( cig==1 )
322         {
323             iread  += ncig;
324             icycle += ncig;
325             continue;
326         }
327         if ( cig==2 )
328         {
329             iref += ncig;
330             continue;
331         }
332         if ( cig==4 )
333         {
334             icycle += ncig;
335             // Soft-clips are present in the sequence, but the position of the read marks a start of non-clipped sequence
336             //   iref += ncig;
337             iread  += ncig;
338             continue;
339         }
340         if ( cig==5 )
341         {
342             icycle += ncig;
343             continue;
344         }
345         if ( cig!=0 )
346             error("TODO: cigar %d, %s\n", cig,bam1_qname(bam_line));
347        
348         if ( ncig+iref > stats->rseq_len )
349             error("FIXME: %d+%d > %d, %s, %s:%d\n",ncig,iref,stats->rseq_len, bam1_qname(bam_line),stats->sam->header->target_name[bam_line->core.tid],bam_line->core.pos+1);
350
351         int im;
352         for (im=0; im<ncig; im++)
353         {
354             uint8_t cread = bam1_seqi(read,iread);
355             uint8_t cref  = stats->rseq_buf[iref];
356
357             // ---------------15
358             // =ACMGRSVTWYHKDBN
359             if ( cread==15 )
360             {
361                 int idx = is_fwd ? icycle : read_len-icycle-1;
362                 if ( idx>stats->max_len )
363                     error("mpc: %d>%d\n",idx,stats->max_len);
364                 idx = idx*stats->nquals;
365                 if ( idx>=stats->nquals*stats->nbases )
366                     error("FIXME: mpc_buf overflow\n");
367                 mpc_buf[idx]++;
368             }
369             else if ( cref && cread && cref!=cread )
370             {
371                 uint8_t qual = quals[iread] + 1;
372                 if ( qual>=stats->nquals )
373                     error("TODO: quality too high %d>=%d\n", quals[iread],stats->nquals);
374
375                 int idx = is_fwd ? icycle : read_len-icycle-1;
376                 if ( idx>stats->max_len )
377                     error("mpc: %d>%d\n",idx,stats->max_len);
378
379                 idx = idx*stats->nquals + qual;
380                 if ( idx>=stats->nquals*stats->nbases )
381                     error("FIXME: mpc_buf overflow\n");
382                 mpc_buf[idx]++;
383             }
384
385             iref++;
386             iread++;
387             icycle++;
388         }
389     }
390 }
391
392 void read_ref_seq(stats_t *stats,int32_t tid,int32_t pos)
393 {
394     khash_t(s) *h;
395     khiter_t iter;
396     faidx1_t val;
397     char *chr, c;
398     faidx_t *fai = stats->fai;
399
400     h = fai->hash;
401     chr = stats->sam->header->target_name[tid];
402
403     // ID of the sequence name
404     iter = kh_get(s, h, chr);
405     if (iter == kh_end(h)) 
406         error("No such reference sequence [%s]?\n", chr);
407     val = kh_value(h, iter);
408
409     // Check the boundaries
410     if (pos >= val.len)
411         error("Was the bam file mapped with the reference sequence supplied?"
412               " A read mapped beyond the end of the chromosome (%s:%d, chromosome length %d).\n", chr,pos,val.len);
413     int size = stats->nref_seq;
414     // The buffer extends beyond the chromosome end. Later the rest will be filled with N's.
415     if (size+pos > val.len) size = val.len-pos;
416
417     // Position the razf reader
418     razf_seek(fai->rz, val.offset + pos / val.line_blen * val.line_len + pos % val.line_blen, SEEK_SET);
419
420     uint8_t *ptr = stats->rseq_buf;
421     int nread = 0;
422     while ( nread<size && razf_read(fai->rz,&c,1) && !fai->rz->z_err )
423     {
424         if ( !isgraph(c) )
425             continue;
426
427         // Conversion between uint8_t coding and ACGT
428         //      -12-4---8-------
429         //      =ACMGRSVTWYHKDBN
430         if ( c=='A' || c=='a' )
431             *ptr = 1;
432         else if ( c=='C' || c=='c' )
433             *ptr = 2;
434         else if ( c=='G' || c=='g' )
435             *ptr = 4;
436         else if ( c=='T' || c=='t' )
437             *ptr = 8;
438         else
439             *ptr = 0;
440         ptr++;
441         nread++;
442     }
443     if ( nread < stats->nref_seq )
444     {
445         memset(ptr,0, stats->nref_seq - nread);
446         nread = stats->nref_seq;
447     }
448     stats->rseq_len = nread;
449     stats->rseq_pos = pos;
450     stats->tid      = tid;
451 }
452
453 float fai_gc_content(stats_t *stats)
454 {
455     uint32_t gc,count,c;
456     int i,size = stats->rseq_len;
457
458     // Count GC content
459     gc = count = 0;
460     for (i=0; i<size; i++)
461     {
462         c = stats->rseq_buf[i];
463         if ( c==2 || c==4 )
464         {
465             gc++;
466             count++;
467         }
468         else if ( c==1 || c==8 )
469             count++;
470     }
471     return count ? (float)gc/count : 0;
472 }
473
474
475 void realloc_buffers(stats_t *stats, int seq_len)
476 {
477     int n = 2*(1 + seq_len - stats->nbases) + stats->nbases;
478
479     stats->quals_1st = realloc(stats->quals_1st, n*stats->nquals*sizeof(uint64_t));
480     if ( !stats->quals_1st )
481         error("Could not realloc buffers, the sequence too long: %d (%ld)\n", seq_len,n*stats->nquals*sizeof(uint64_t));
482     memset(stats->quals_1st + stats->nbases*stats->nquals, 0, (n-stats->nbases)*stats->nquals*sizeof(uint64_t));
483
484     stats->quals_2nd = realloc(stats->quals_2nd, n*stats->nquals*sizeof(uint64_t));
485     if ( !stats->quals_2nd )
486         error("Could not realloc buffers, the sequence too long: %d (2x%ld)\n", seq_len,n*stats->nquals*sizeof(uint64_t));
487     memset(stats->quals_2nd + stats->nbases*stats->nquals, 0, (n-stats->nbases)*stats->nquals*sizeof(uint64_t));
488
489     if ( stats->mpc_buf )
490     {
491         stats->mpc_buf = realloc(stats->mpc_buf, n*stats->nquals*sizeof(uint64_t));
492         if ( !stats->mpc_buf )
493             error("Could not realloc buffers, the sequence too long: %d (%ld)\n", seq_len,n*stats->nquals*sizeof(uint64_t));
494         memset(stats->mpc_buf + stats->nbases*stats->nquals, 0, (n-stats->nbases)*stats->nquals*sizeof(uint64_t));
495     }
496
497     stats->acgt_cycles = realloc(stats->acgt_cycles, n*4*sizeof(uint64_t));
498     if ( !stats->acgt_cycles )
499         error("Could not realloc buffers, the sequence too long: %d (%ld)\n", seq_len,n*4*sizeof(uint64_t));
500     memset(stats->acgt_cycles + stats->nbases*4, 0, (n-stats->nbases)*4*sizeof(uint64_t));
501
502     stats->read_lengths = realloc(stats->read_lengths, n*sizeof(uint64_t));
503     if ( !stats->read_lengths )
504         error("Could not realloc buffers, the sequence too long: %d (%ld)\n", seq_len,n*sizeof(uint64_t));
505     memset(stats->read_lengths + stats->nbases, 0, (n-stats->nbases)*sizeof(uint64_t));
506
507     stats->insertions = realloc(stats->insertions, n*sizeof(uint64_t));
508     if ( !stats->insertions )
509         error("Could not realloc buffers, the sequence too long: %d (%ld)\n", seq_len,n*sizeof(uint64_t));
510     memset(stats->insertions + stats->nbases, 0, (n-stats->nbases)*sizeof(uint64_t));
511
512     stats->deletions = realloc(stats->deletions, n*sizeof(uint64_t));
513     if ( !stats->deletions )
514         error("Could not realloc buffers, the sequence too long: %d (%ld)\n", seq_len,n*sizeof(uint64_t));
515     memset(stats->deletions + stats->nbases, 0, (n-stats->nbases)*sizeof(uint64_t));
516
517     stats->ins_cycles = realloc(stats->ins_cycles, n*sizeof(uint64_t));
518     if ( !stats->ins_cycles )
519         error("Could not realloc buffers, the sequence too long: %d (%ld)\n", seq_len,n*sizeof(uint64_t));
520     memset(stats->ins_cycles + stats->nbases, 0, (n-stats->nbases)*sizeof(uint64_t));
521
522     stats->del_cycles = realloc(stats->del_cycles, n*sizeof(uint64_t));
523     if ( !stats->del_cycles )
524         error("Could not realloc buffers, the sequence too long: %d (%ld)\n", seq_len,n*sizeof(uint64_t));
525     memset(stats->del_cycles + stats->nbases, 0, (n-stats->nbases)*sizeof(uint64_t));
526
527     stats->nbases = n;
528
529     // Realloc the coverage distribution buffer
530     int *rbuffer = calloc(sizeof(int),seq_len*5);
531     n = stats->cov_rbuf.size-stats->cov_rbuf.start;
532     memcpy(rbuffer,stats->cov_rbuf.buffer+stats->cov_rbuf.start,n);
533     if ( stats->cov_rbuf.start>1 )
534         memcpy(rbuffer+n,stats->cov_rbuf.buffer,stats->cov_rbuf.start);
535     stats->cov_rbuf.start = 0;
536     free(stats->cov_rbuf.buffer);
537     stats->cov_rbuf.buffer = rbuffer;
538     stats->cov_rbuf.size = seq_len*5;
539 }
540
541 void collect_stats(bam1_t *bam_line, stats_t *stats)
542 {
543     if ( stats->rmdup && IS_DUP(bam_line) )
544         return;
545
546     int seq_len = bam_line->core.l_qseq;
547     if ( !seq_len ) return;
548     if ( stats->filter_readlen!=-1 && seq_len!=stats->filter_readlen ) return;
549     if ( seq_len >= stats->nbases )
550         realloc_buffers(stats,seq_len);
551     if ( stats->max_len<seq_len )
552         stats->max_len = seq_len;
553
554     stats->read_lengths[seq_len]++;
555
556     // Count GC and ACGT per cycle
557     uint8_t base, *seq  = bam1_seq(bam_line);
558     int gc_count  = 0;
559     int i;
560     int reverse = IS_REVERSE(bam_line);
561     for (i=0; i<seq_len; i++)
562     {
563         // Conversion from uint8_t coding to ACGT
564         //      -12-4---8-------
565         //      =ACMGRSVTWYHKDBN
566         //       01 2   3
567         base = bam1_seqi(seq,i);
568         base /= 2;
569         if ( base==1 || base==2 ) gc_count++;
570         else if ( base>2 ) base=3;
571         if ( 4*(reverse ? seq_len-i-1 : i) + base >= stats->nbases*4 ) 
572             error("FIXME: acgt_cycles\n");
573         stats->acgt_cycles[ 4*(reverse ? seq_len-i-1 : i) + base ]++;
574     }
575     int gc_idx_min = gc_count*(stats->ngc-1)/seq_len;
576     int gc_idx_max = (gc_count+1)*(stats->ngc-1)/seq_len;
577     if ( gc_idx_max >= stats->ngc ) gc_idx_max = stats->ngc - 1;
578
579     // Determine which array (1st or 2nd read) will these stats go to,
580     //  trim low quality bases from end the same way BWA does, 
581     //  fill GC histogram
582     uint64_t *quals;
583     uint8_t *bam_quals = bam1_qual(bam_line);
584     if ( bam_line->core.flag&BAM_FREAD2 )
585     {
586         quals  = stats->quals_2nd;
587         stats->nreads_2nd++;
588         for (i=gc_idx_min; i<gc_idx_max; i++)
589             stats->gc_2nd[i]++;
590     }
591     else
592     {
593         quals = stats->quals_1st;
594         stats->nreads_1st++;
595         for (i=gc_idx_min; i<gc_idx_max; i++)
596             stats->gc_1st[i]++;
597     }
598     if ( stats->trim_qual>0 ) 
599         stats->nbases_trimmed += bwa_trim_read(stats->trim_qual, bam_quals, seq_len, reverse);
600
601     // Quality histogram and average quality
602     for (i=0; i<seq_len; i++)
603     {
604         uint8_t qual = bam_quals[ reverse ? seq_len-i-1 : i];
605         if ( qual>=stats->nquals )
606             error("TODO: quality too high %d>=%d\n", quals[i],stats->nquals);
607         if ( qual>stats->max_qual )
608             stats->max_qual = qual;
609
610         quals[ i*stats->nquals+qual ]++;
611         stats->sum_qual += qual;
612     }
613
614     // Look at the flags and increment appropriate counters (mapped, paired, etc)
615     if ( IS_UNMAPPED(bam_line) )
616         stats->nreads_unmapped++;
617     else
618     {
619         if ( !bam_line->core.qual )
620             stats->nreads_mq0++;
621
622         count_indels(stats,bam_line);
623
624         // The insert size is tricky, because for long inserts the libraries are
625         // prepared differently and the pairs point in other direction. BWA does
626         // not set the paired flag for them.  Similar thing is true also for 454
627         // reads. Therefore, do the insert size stats for all mapped reads.
628         int32_t isize = bam_line->core.isize;
629         if ( isize<0 ) isize = -isize;
630         if ( IS_PAIRED(bam_line) && isize!=0 )
631         {
632             stats->nreads_paired++;
633             if ( isize >= stats->nisize )
634                 isize=stats->nisize-1;
635
636             int pos_fst = bam_line->core.mpos - bam_line->core.pos;
637             int is_fst  = IS_READ1(bam_line) ? 1 : -1;
638             int is_fwd  = IS_REVERSE(bam_line) ? -1 : 1;
639             int is_mfwd = IS_MATE_REVERSE(bam_line) ? -1 : 1;
640
641             if ( is_fwd*is_mfwd>0 )
642                 stats->isize_other[isize]++;
643             else if ( is_fst*pos_fst>0 )
644             {
645                 if ( is_fst*is_fwd>0 )
646                     stats->isize_inward[isize]++;
647                 else
648                     stats->isize_outward[isize]++;
649             }
650             else if ( is_fst*pos_fst<0 )
651             {
652                 if ( is_fst*is_fwd>0 )
653                     stats->isize_outward[isize]++;
654                 else
655                     stats->isize_inward[isize]++;
656             }
657         }
658         else
659             stats->nreads_unpaired++;
660
661         // Number of mismatches 
662         uint8_t *nm = bam_aux_get(bam_line,"NM");
663         if (nm) 
664             stats->nmismatches += bam_aux2i(nm);
665
666         // Number of mapped bases from cigar 
667         if ( bam_line->core.n_cigar == 0) 
668             error("FIXME: mapped read with no cigar?\n");
669         int readlen = seq_len;
670         for (i=0; i<bam_line->core.n_cigar; i++) 
671         {
672             // Conversion from uint32_t to MIDNSHP
673             //  01--4--
674             //  MIDNSHP
675             if ( (bam1_cigar(bam_line)[i]&BAM_CIGAR_MASK)==0 || (bam1_cigar(bam_line)[i]&BAM_CIGAR_MASK)==1 )
676                 stats->nbases_mapped_cigar += bam1_cigar(bam_line)[i]>>BAM_CIGAR_SHIFT;
677
678             if ( (bam1_cigar(bam_line)[i]&BAM_CIGAR_MASK)==2 )
679                 readlen += bam1_cigar(bam_line)[i]>>BAM_CIGAR_SHIFT;
680         }
681         stats->nbases_mapped += seq_len;
682
683         if ( stats->tid==bam_line->core.tid && bam_line->core.pos<stats->pos )
684             stats->is_sorted = 0;
685         stats->pos = bam_line->core.pos;
686
687         if ( stats->is_sorted )
688         {
689             if ( stats->tid==-1 || stats->tid!=bam_line->core.tid )
690                 round_buffer_flush(stats,-1);
691
692             // Mismatches per cycle and GC-depth graph
693             if ( stats->fai )
694             {
695                 // First pass, new chromosome or sequence spanning beyond the end of the buffer 
696                 if ( stats->rseq_pos==-1 || stats->tid != bam_line->core.tid || stats->rseq_pos+stats->rseq_len < bam_line->core.pos+readlen )
697                 {
698                     read_ref_seq(stats,bam_line->core.tid,bam_line->core.pos);
699
700                     // Get the reference GC content for this bin. Note that in the current implementation
701                     //  the GCD bins overlap by seq_len-1. Because the bin size is by default 20k and the 
702                     //  expected read length only ~100bp, it shouldn't really matter.
703                     stats->gcd_pos = bam_line->core.pos;
704                     stats->igcd++;
705
706                     if ( stats->igcd >= stats->ngcd )
707                         error("The genome too long or the GCD bin overlaps too big [%ud]\n", stats->igcd);
708
709                     stats->gcd[ stats->igcd ].gc = fai_gc_content(stats);
710                 }
711                 count_mismatches_per_cycle(stats,bam_line);
712             }
713             else if ( stats->gcd_pos==-1 || stats->tid != bam_line->core.tid || bam_line->core.pos - stats->gcd_pos > stats->gcd_bin_size )
714             {
715                 // First pass or a new chromosome
716                 stats->tid     = bam_line->core.tid;
717                 stats->gcd_pos = bam_line->core.pos;
718                 stats->igcd++;
719
720                 if ( stats->igcd >= stats->ngcd )
721                 {
722                     uint32_t n = 2*(1 + stats->ngcd);
723                     stats->gcd = realloc(stats->gcd, n*sizeof(gc_depth_t));
724                     if ( !stats->gcd )
725                         error("Could not realloc GCD buffer, too many chromosomes or the genome too long?? [%u %u]\n", stats->ngcd,n);
726                     memset(&(stats->gcd[stats->ngcd]),0,(n-stats->ngcd)*sizeof(gc_depth_t));
727                 }
728             }
729
730             stats->gcd[ stats->igcd ].depth++;
731             // When no reference sequence is given, approximate the GC from the read (much shorter window, but otherwise OK)
732             if ( !stats->fai )
733                 stats->gcd[ stats->igcd ].gc += (float) gc_count / seq_len;
734
735             // Coverage distribution graph
736             round_buffer_flush(stats,bam_line->core.pos);
737             round_buffer_insert_read(&(stats->cov_rbuf),bam_line->core.pos,bam_line->core.pos+seq_len-1);
738         }
739     }
740
741     stats->total_len += seq_len;
742     if ( IS_DUP(bam_line) )
743     {
744         stats->total_len_dup += seq_len;
745         stats->nreads_dup++;
746     }
747 }
748
749 // Sort by GC and depth
750 #define GCD_t(x) ((gc_depth_t *)x)
751 static int gcd_cmp(const void *a, const void *b)
752 {
753     if ( GCD_t(a)->gc < GCD_t(b)->gc ) return -1;
754     if ( GCD_t(a)->gc > GCD_t(b)->gc ) return 1;
755     if ( GCD_t(a)->depth < GCD_t(b)->depth ) return -1;
756     if ( GCD_t(a)->depth > GCD_t(b)->depth ) return 1;
757     return 0;
758 }
759 #undef GCD_t
760
761 float gcd_percentile(gc_depth_t *gcd, int N, int p)
762 {
763     float n,d;
764     int k;
765
766     n = p*(N+1)/100;
767     k = n;
768     if ( k<=0 ) 
769         return gcd[0].depth;
770     if ( k>=N ) 
771         return gcd[N-1].depth;
772
773     d = n - k;
774     return gcd[k-1].depth + d*(gcd[k].depth - gcd[k-1].depth);
775 }
776
777 void output_stats(stats_t *stats)
778 {
779     // Calculate average insert size and standard deviation (from the main bulk data only)
780     int isize, ibulk=0;
781     uint64_t nisize=0, nisize_inward=0, nisize_outward=0, nisize_other=0;
782     for (isize=1; isize<stats->nisize; isize++)
783     {
784         // Each pair was counted twice
785         stats->isize_inward[isize] *= 0.5;
786         stats->isize_outward[isize] *= 0.5;
787         stats->isize_other[isize] *= 0.5;
788
789         nisize_inward += stats->isize_inward[isize];
790         nisize_outward += stats->isize_outward[isize];
791         nisize_other += stats->isize_other[isize];
792         nisize += stats->isize_inward[isize] + stats->isize_outward[isize] + stats->isize_other[isize];
793     }
794
795     double bulk=0, avg_isize=0, sd_isize=0;
796     for (isize=1; isize<stats->nisize; isize++)
797     {
798         bulk += stats->isize_inward[isize] + stats->isize_outward[isize] + stats->isize_other[isize];
799         avg_isize += isize * (stats->isize_inward[isize] + stats->isize_outward[isize] + stats->isize_other[isize]);
800
801         if ( bulk/nisize > stats->isize_main_bulk )
802         {
803             ibulk  = isize+1;
804             nisize = bulk;
805             break;
806         }
807     }
808     avg_isize /= nisize ? nisize : 1;
809     for (isize=1; isize<ibulk; isize++)
810         sd_isize += (stats->isize_inward[isize] + stats->isize_outward[isize] + stats->isize_other[isize]) * (isize-avg_isize)*(isize-avg_isize) / nisize;
811     sd_isize = sqrt(sd_isize);
812
813
814     printf("# This file was produced by bamcheck (%s)\n",BAMCHECK_VERSION);
815     printf("# The command line was:  %s",stats->argv[0]);
816     int i;
817     for (i=1; i<stats->argc; i++)
818         printf(" %s",stats->argv[i]);
819     printf("\n");
820     printf("# Summary Numbers. Use `grep ^SN | cut -f 2-` to extract this part.\n");
821     printf("SN\tsequences:\t%ld\n", (long)(stats->nreads_1st+stats->nreads_2nd));
822     printf("SN\tis paired:\t%d\n", stats->nreads_1st&&stats->nreads_2nd ? 1 : 0);
823     printf("SN\tis sorted:\t%d\n", stats->is_sorted ? 1 : 0);
824     printf("SN\t1st fragments:\t%ld\n", (long)stats->nreads_1st);
825     printf("SN\tlast fragments:\t%ld\n", (long)stats->nreads_2nd);
826     printf("SN\treads mapped:\t%ld\n", (long)(stats->nreads_paired+stats->nreads_unpaired));
827     printf("SN\treads unmapped:\t%ld\n", (long)stats->nreads_unmapped);
828     printf("SN\treads unpaired:\t%ld\n", (long)stats->nreads_unpaired);
829     printf("SN\treads paired:\t%ld\n", (long)stats->nreads_paired);
830     printf("SN\treads duplicated:\t%ld\n", (long)stats->nreads_dup);
831     printf("SN\treads MQ0:\t%ld\n", (long)stats->nreads_mq0);
832     printf("SN\ttotal length:\t%ld\n", (long)stats->total_len);
833     printf("SN\tbases mapped:\t%ld\n", (long)stats->nbases_mapped);
834     printf("SN\tbases mapped (cigar):\t%ld\n", (long)stats->nbases_mapped_cigar);
835     printf("SN\tbases trimmed:\t%ld\n", (long)stats->nbases_trimmed);
836     printf("SN\tbases duplicated:\t%ld\n", (long)stats->total_len_dup);
837     printf("SN\tmismatches:\t%ld\n", (long)stats->nmismatches);
838     printf("SN\terror rate:\t%e\n", (float)stats->nmismatches/stats->nbases_mapped_cigar);
839     float avg_read_length = (stats->nreads_1st+stats->nreads_2nd)?stats->total_len/(stats->nreads_1st+stats->nreads_2nd):0;
840     printf("SN\taverage length:\t%.0f\n", avg_read_length);
841     printf("SN\tmaximum length:\t%d\n", stats->max_len);
842     printf("SN\taverage quality:\t%.1f\n", stats->total_len?stats->sum_qual/stats->total_len:0);
843     printf("SN\tinsert size average:\t%.1f\n", avg_isize);
844     printf("SN\tinsert size standard deviation:\t%.1f\n", sd_isize);
845     printf("SN\tinward oriented pairs:\t%ld\n", (long)nisize_inward);
846     printf("SN\toutward oriented pairs:\t%ld\n", (long)nisize_outward);
847     printf("SN\tpairs with other orientation:\t%ld\n", (long)nisize_other);
848
849     int ibase,iqual;
850     if ( stats->max_len<stats->nbases ) stats->max_len++;
851     if ( stats->max_qual+1<stats->nquals ) stats->max_qual++;
852     printf("# First Fragment Qualitites. Use `grep ^FFQ | cut -f 2-` to extract this part.\n");
853     printf("# Columns correspond to qualities and rows to cycles. First column is the cycle number.\n");
854     for (ibase=0; ibase<stats->max_len; ibase++)
855     {
856         printf("FFQ\t%d",ibase+1);
857         for (iqual=0; iqual<=stats->max_qual; iqual++)
858         {
859             printf("\t%ld", (long)stats->quals_1st[ibase*stats->nquals+iqual]);
860         }
861         printf("\n");
862     }
863     printf("# Last Fragment Qualitites. Use `grep ^LFQ | cut -f 2-` to extract this part.\n");
864     printf("# Columns correspond to qualities and rows to cycles. First column is the cycle number.\n");
865     for (ibase=0; ibase<stats->max_len; ibase++)
866     {
867         printf("LFQ\t%d",ibase+1);
868         for (iqual=0; iqual<=stats->max_qual; iqual++)
869         {
870             printf("\t%ld", (long)stats->quals_2nd[ibase*stats->nquals+iqual]);
871         }
872         printf("\n");
873     }
874     if ( stats->mpc_buf )
875     {
876         printf("# Mismatches per cycle and quality. Use `grep ^MPC | cut -f 2-` to extract this part.\n");
877         printf("# Columns correspond to qualities, rows to cycles. First column is the cycle number, second\n");
878         printf("# is the number of N's and the rest is the number of mismatches\n");
879         for (ibase=0; ibase<stats->max_len; ibase++)
880         {
881             printf("MPC\t%d",ibase+1);
882             for (iqual=0; iqual<=stats->max_qual; iqual++)
883             {
884                 printf("\t%ld", (long)stats->mpc_buf[ibase*stats->nquals+iqual]);
885             }
886             printf("\n");
887         }
888     }
889     printf("# GC Content of first fragments. Use `grep ^GCF | cut -f 2-` to extract this part.\n");
890     int ibase_prev = 0;
891     for (ibase=0; ibase<stats->ngc; ibase++)
892     {
893         if ( stats->gc_1st[ibase]==stats->gc_1st[ibase_prev] ) continue;
894         printf("GCF\t%.2f\t%ld\n", (ibase+ibase_prev)*0.5*100./(stats->ngc-1), (long)stats->gc_1st[ibase_prev]);
895         ibase_prev = ibase;
896     }
897     printf("# GC Content of last fragments. Use `grep ^GCL | cut -f 2-` to extract this part.\n");
898     ibase_prev = 0;
899     for (ibase=0; ibase<stats->ngc; ibase++)
900     {
901         if ( stats->gc_2nd[ibase]==stats->gc_2nd[ibase_prev] ) continue;
902         printf("GCL\t%.2f\t%ld\n", (ibase+ibase_prev)*0.5*100./(stats->ngc-1), (long)stats->gc_2nd[ibase_prev]);
903         ibase_prev = ibase;
904     }
905     printf("# ACGT content per cycle. Use `grep ^GCC | cut -f 2-` to extract this part. The columns are: cycle, and A,C,G,T counts [%%]\n");
906     for (ibase=0; ibase<stats->max_len; ibase++)
907     {
908         uint64_t *ptr = &(stats->acgt_cycles[ibase*4]);
909         uint64_t  sum = ptr[0]+ptr[1]+ptr[2]+ptr[3];
910         if ( ! sum ) continue;
911         printf("GCC\t%d\t%.2f\t%.2f\t%.2f\t%.2f\n", ibase,100.*ptr[0]/sum,100.*ptr[1]/sum,100.*ptr[2]/sum,100.*ptr[3]/sum);
912     }
913     printf("# Insert sizes. Use `grep ^IS | cut -f 2-` to extract this part. The columns are: pairs total, inward oriented pairs, outward oriented pairs, other pairs\n");
914     for (isize=1; isize<ibulk; isize++)
915         printf("IS\t%d\t%ld\t%ld\t%ld\t%ld\n", isize, (long)(stats->isize_inward[isize]+stats->isize_outward[isize]+stats->isize_other[isize]),
916             (long)stats->isize_inward[isize], (long)stats->isize_outward[isize], (long)stats->isize_other[isize]);
917
918     printf("# Read lengths. Use `grep ^RL | cut -f 2-` to extract this part. The columns are: read length, count\n");
919     int ilen;
920     for (ilen=0; ilen<stats->max_len; ilen++)
921     {
922         if ( stats->read_lengths[ilen]>0 )
923             printf("RL\t%d\t%ld\n", ilen, (long)stats->read_lengths[ilen]);
924     }
925
926     printf("# Indel distribution. Use `grep ^ID | cut -f 2-` to extract this part. The columns are: length, number of insertions, number of deletions\n");
927     for (ilen=0; ilen<stats->nindels; ilen++)
928     {
929         if ( stats->insertions[ilen]>0 || stats->deletions[ilen]>0 )
930             printf("ID\t%d\t%ld\t%ld\n", ilen+1, (long)stats->insertions[ilen], (long)stats->deletions[ilen]);
931     }
932
933     printf("# Indels per cycle. Use `grep ^IC | cut -f 2-` to extract this part. The columns are: cycle, number of insertions, number of deletions\n");
934     for (ilen=0; ilen<stats->nbases; ilen++)
935     {
936         if ( stats->ins_cycles[ilen]>0 || stats->del_cycles[ilen]>0 )
937             printf("IC\t%d\t%ld\t%ld\n", ilen+1, (long)stats->ins_cycles[ilen], (long)stats->del_cycles[ilen]);
938     }
939
940     printf("# Coverage distribution. Use `grep ^COV | cut -f 2-` to extract this part.\n");
941     printf("COV\t[<%d]\t%d\t%ld\n",stats->cov_min,stats->cov_min-1, (long)stats->cov[0]);
942     int icov;
943     for (icov=1; icov<stats->ncov-1; icov++)
944         printf("COV\t[%d-%d]\t%d\t%ld\n",stats->cov_min + (icov-1)*stats->cov_step, stats->cov_min + icov*stats->cov_step-1,stats->cov_min + icov*stats->cov_step-1, (long)stats->cov[icov]);
945     printf("COV\t[%d<]\t%d\t%ld\n",stats->cov_min + (stats->ncov-2)*stats->cov_step-1,stats->cov_min + (stats->ncov-2)*stats->cov_step-1, (long)stats->cov[stats->ncov-1]);
946
947
948     // Calculate average GC content, then sort by GC and depth
949     printf("# GC-depth. Use `grep ^GCD | cut -f 2-` to extract this part. The columns are: GC%%, unique sequence percentiles, 10th, 25th, 50th, 75th and 90th depth percentile\n");
950     uint32_t igcd;
951     for (igcd=0; igcd<stats->igcd; igcd++)
952     {
953         if ( stats->fai )
954             stats->gcd[igcd].gc = round(100. * stats->gcd[igcd].gc);
955         else
956             if ( stats->gcd[igcd].depth ) 
957                 stats->gcd[igcd].gc = round(100. * stats->gcd[igcd].gc / stats->gcd[igcd].depth);
958     }
959     qsort(stats->gcd, stats->igcd+1, sizeof(gc_depth_t), gcd_cmp);
960     igcd = 0;
961     while ( igcd < stats->igcd )
962     {
963         // Calculate percentiles (10,25,50,75,90th) for the current GC content and print
964         uint32_t nbins=0, itmp=igcd;
965         float gc = stats->gcd[igcd].gc;
966         while ( itmp<stats->igcd && fabs(stats->gcd[itmp].gc-gc)<0.1 )
967         {
968             nbins++;
969             itmp++;
970         }
971         printf("GCD\t%.1f\t%.3f\t%.3f\t%.3f\t%.3f\t%.3f\t%.3f\n", gc, (igcd+nbins+1)*100./(stats->igcd+1),
972                 gcd_percentile(&(stats->gcd[igcd]),nbins,10) *avg_read_length/stats->gcd_bin_size,
973                 gcd_percentile(&(stats->gcd[igcd]),nbins,25) *avg_read_length/stats->gcd_bin_size, 
974                 gcd_percentile(&(stats->gcd[igcd]),nbins,50) *avg_read_length/stats->gcd_bin_size, 
975                 gcd_percentile(&(stats->gcd[igcd]),nbins,75) *avg_read_length/stats->gcd_bin_size, 
976                 gcd_percentile(&(stats->gcd[igcd]),nbins,90) *avg_read_length/stats->gcd_bin_size 
977               );
978         igcd += nbins;
979     }
980 }
981
982 void bam_init_header_hash(bam_header_t *header);
983
984 size_t getline(char **line, size_t *n, FILE *fp)
985 {
986     if (line == NULL || n == NULL || fp == NULL)
987     {
988         errno = EINVAL;
989         return -1;
990     }
991     if (*n==0 || !*line)
992     {
993         *line = NULL;
994         *n = 0;
995     }
996
997     size_t nread=0;
998     int c;
999     while ((c=getc(fp))!= EOF && c!='\n')
1000     {
1001         if ( ++nread>=*n )
1002         {
1003             *n += 255;
1004             *line = realloc(*line, sizeof(char)*(*n));
1005         }
1006         (*line)[nread-1] = c;
1007     }
1008     if ( nread>=*n )
1009     {
1010         *n += 255;
1011         *line = realloc(*line, sizeof(char)*(*n));
1012     }
1013     (*line)[nread] = 0;
1014     return nread>0 ? nread : -1;
1015
1016 }
1017
1018 void init_regions(stats_t *stats, char *file)
1019 {
1020     khiter_t iter;
1021     khash_t(str) *header_hash;
1022
1023     bam_init_header_hash(stats->sam->header);
1024     header_hash = (khash_t(str)*)stats->sam->header->hash;
1025
1026     FILE *fp = fopen(file,"r");
1027     if ( !fp ) error("%s: %s\n",file,strerror(errno));
1028
1029     char *line = NULL;
1030     size_t len = 0;
1031     ssize_t nread;
1032     int warned = 0;
1033     int prev_tid=-1, prev_pos=-1;
1034     while ((nread = getline(&line, &len, fp)) != -1) 
1035     {
1036         if ( line[0] == '#' ) continue;
1037
1038         int i = 0;
1039         while ( i<nread && !isspace(line[i]) ) i++;
1040         if ( i>=nread ) error("Could not parse the file: %s [%s]\n", file,line);
1041         line[i] = 0;
1042
1043         iter = kh_get(str, header_hash, line);
1044         int tid = kh_val(header_hash, iter);
1045         if ( iter == kh_end(header_hash) )
1046         {
1047             if ( !warned )
1048                 fprintf(stderr,"Warning: Some sequences not present in the BAM (%s)\n", line);
1049             warned = 1;
1050             continue;
1051         }
1052
1053         if ( tid >= stats->nregions )
1054         {
1055             stats->regions = realloc(stats->regions,sizeof(regions_t)*(stats->nregions+100));
1056             int j;
1057             for (j=stats->nregions; j<stats->nregions+100; j++)
1058             {
1059                 stats->regions[j].npos = stats->regions[j].mpos = stats->regions[j].cpos = 0;
1060                 stats->regions[j].pos = NULL;
1061             }
1062             stats->nregions += 100;
1063         }
1064         int npos = stats->regions[tid].npos;
1065         if ( npos >= stats->regions[tid].mpos )
1066         {
1067             stats->regions[tid].mpos += 1000;
1068             stats->regions[tid].pos = realloc(stats->regions[tid].pos,sizeof(pos_t)*stats->regions[tid].mpos);
1069         }
1070
1071         if ( (sscanf(line+i+1,"%d %d",&stats->regions[tid].pos[npos].from,&stats->regions[tid].pos[npos].to))!=2 ) error("Could not parse the region [%s]\n");
1072         if ( prev_tid==-1 || prev_tid!=tid )
1073         {
1074             prev_tid = tid;
1075             prev_pos = stats->regions[tid].pos[npos].from;
1076         }
1077         if ( prev_pos>stats->regions[tid].pos[npos].from )
1078             error("The positions are not in chromosomal order (%s:%d comes after %d)\n", line,stats->regions[tid].pos[npos].from,prev_pos);
1079         stats->regions[tid].npos++;
1080     }
1081     if (line) free(line);
1082     fclose(fp);
1083 }
1084
1085 void destroy_regions(stats_t *stats)
1086 {
1087     int i;
1088     for (i=0; i<stats->nregions; i++)
1089     {
1090         if ( !stats->regions[i].mpos ) continue;
1091         free(stats->regions[i].pos);
1092     }
1093     if ( stats->regions ) free(stats->regions);
1094 }
1095
1096 static int fetch_read(const bam1_t *bam_line, void *data)
1097 {
1098     collect_stats((bam1_t*)bam_line,(stats_t*)data);
1099     return 1;
1100 }
1101
1102
1103 void error(const char *format, ...)
1104 {
1105     if ( !format )
1106     {
1107         printf("Version: %s\n", BAMCHECK_VERSION);
1108         printf("About: The program collects statistics from BAM files. The output can be visualized using plot-bamcheck.\n");
1109         printf("Usage: bamcheck [OPTIONS] file.bam\n");
1110         printf("       bamcheck [OPTIONS] file.bam chr:from-to\n");
1111         printf("Options:\n");
1112         printf("    -c, --coverage <int>,<int>,<int>    Coverage distribution min,max,step [1,1000,1]\n");
1113         printf("    -d, --remove-dups                   Exlude from statistics reads marked as duplicates\n");
1114         printf("    -h, --help                          This help message\n");
1115         printf("    -i, --insert-size <int>             Maximum insert size [8000]\n");
1116         printf("    -l, --read-length <int>             Include in the statistics only reads with the given read length []\n");
1117         printf("    -m, --most-inserts <float>          Report only the main part of inserts [0.99]\n");
1118         printf("    -q, --trim-quality <int>            The BWA trimming parameter [0]\n");
1119         printf("    -r, --ref-seq <file>                Reference sequence (required for GC-depth calculation).\n");
1120         printf("    -t, --target-regions <file>         Do stats in these regions only. Tab-delimited file chr,from,to, 1-based, inclusive.\n");
1121         printf("    -s, --sam                           Input is SAM\n");
1122         printf("\n");
1123     }
1124     else
1125     {
1126         va_list ap;
1127         va_start(ap, format);
1128         vfprintf(stderr, format, ap);
1129         va_end(ap);
1130     }
1131     exit(-1);
1132 }
1133
1134 int main(int argc, char *argv[])
1135 {
1136     char *targets = NULL;
1137     char *bam_fname = NULL;
1138     samfile_t *sam = NULL;
1139     char in_mode[5];
1140
1141     stats_t *stats = calloc(1,sizeof(stats_t));
1142     stats->ngc    = 200;
1143     stats->nquals = 95;
1144     stats->nbases = 300;
1145     stats->nisize = 8000;
1146     stats->max_len   = 30;
1147     stats->max_qual  = 40;
1148     stats->isize_main_bulk = 0.99;   // There are always outliers at the far end
1149     stats->gcd_bin_size = 20000;
1150     stats->ngcd         = 3e5;     // 300k of 20k bins is enough to hold a genome 6Gbp big
1151     stats->nref_seq     = stats->gcd_bin_size;
1152     stats->rseq_pos     = -1;
1153     stats->tid = stats->gcd_pos = -1;
1154     stats->is_sorted = 1;
1155     stats->cov_min  = 1;
1156     stats->cov_max  = 1000;
1157     stats->cov_step = 1;
1158     stats->argc = argc;
1159     stats->argv = argv;
1160     stats->filter_readlen = -1;
1161     stats->nindels = stats->nbases;
1162
1163     strcpy(in_mode, "rb");
1164
1165     static struct option loptions[] = 
1166     {
1167         {"help",0,0,'h'},
1168         {"remove-dups",0,0,'d'},
1169         {"sam",0,0,'s'},
1170         {"ref-seq",0,0,'r'},
1171         {"coverage",0,0,'c'},
1172         {"read-length",0,0,'l'},
1173         {"insert-size",0,0,'i'},
1174         {"most-inserts",0,0,'m'},
1175         {"trim-quality",0,0,'q'},
1176         {"target-regions",0,0,'t'},
1177         {0,0,0,0}
1178     };
1179     int opt;
1180     while ( (opt=getopt_long(argc,argv,"?hdsr:c:l:i:t:m:q:",loptions,NULL))>0 )
1181     {
1182         switch (opt)
1183         {
1184             case 'd': stats->rmdup=1; break;
1185             case 's': strcpy(in_mode, "r"); break;
1186             case 'r': stats->fai = fai_load(optarg); 
1187                       if (stats->fai==0) 
1188                           error("Could not load faidx: %s\n", optarg); 
1189                       break;
1190             case 'c': if ( sscanf(optarg,"%d,%d,%d",&stats->cov_min,&stats->cov_max,&stats->cov_step)!= 3 ) 
1191                           error("Unable to parse -c %s\n", optarg); 
1192                       break;
1193             case 'l': stats->filter_readlen = atoi(optarg); break;
1194             case 'i': stats->nisize = atoi(optarg); break;
1195             case 'm': stats->isize_main_bulk = atof(optarg); break;
1196             case 'q': stats->trim_qual = atoi(optarg); break;
1197             case 't': targets = optarg; break;
1198             case '?': 
1199             case 'h': error(NULL);
1200             default: error("Unknown argument: %s\n", optarg);
1201         }
1202     }
1203     if ( optind<argc )
1204         bam_fname = argv[optind++];
1205
1206     if ( !bam_fname )
1207     {
1208         if ( isatty(fileno((FILE *)stdin)) )
1209             error(NULL);
1210         bam_fname = "-";
1211     }
1212
1213     // Init structures
1214     //  .. coverage bins and round buffer
1215     if ( stats->cov_step > stats->cov_max - stats->cov_min + 1 )
1216     {
1217         stats->cov_step = stats->cov_max - stats->cov_min;
1218         if ( stats->cov_step <= 0 )
1219             stats->cov_step = 1;
1220     }
1221     stats->ncov = 3 + (stats->cov_max-stats->cov_min) / stats->cov_step;
1222     stats->cov_max = stats->cov_min + ((stats->cov_max-stats->cov_min)/stats->cov_step +1)*stats->cov_step - 1;
1223     stats->cov = calloc(sizeof(uint64_t),stats->ncov);
1224     stats->cov_rbuf.size = stats->nbases*5;
1225     stats->cov_rbuf.buffer = calloc(sizeof(int32_t),stats->cov_rbuf.size);
1226     // .. bam
1227     if ((sam = samopen(bam_fname, in_mode, NULL)) == 0) 
1228         error("Failed to open: %s\n", bam_fname);
1229     stats->sam = sam;
1230     bam1_t *bam_line = bam_init1();
1231     // .. arrays
1232     stats->quals_1st     = calloc(stats->nquals*stats->nbases,sizeof(uint64_t));
1233     stats->quals_2nd     = calloc(stats->nquals*stats->nbases,sizeof(uint64_t));
1234     stats->gc_1st        = calloc(stats->ngc,sizeof(uint64_t));
1235     stats->gc_2nd        = calloc(stats->ngc,sizeof(uint64_t));
1236     stats->isize_inward  = calloc(stats->nisize,sizeof(uint64_t));
1237     stats->isize_outward = calloc(stats->nisize,sizeof(uint64_t));
1238     stats->isize_other   = calloc(stats->nisize,sizeof(uint64_t));
1239     stats->gcd           = calloc(stats->ngcd,sizeof(gc_depth_t));
1240     stats->rseq_buf      = calloc(stats->nref_seq,sizeof(uint8_t));
1241     stats->mpc_buf       = stats->fai ? calloc(stats->nquals*stats->nbases,sizeof(uint64_t)) : NULL;
1242     stats->acgt_cycles   = calloc(4*stats->nbases,sizeof(uint64_t));
1243     stats->read_lengths  = calloc(stats->nbases,sizeof(uint64_t));
1244     stats->insertions    = calloc(stats->nbases,sizeof(uint64_t));
1245     stats->deletions     = calloc(stats->nbases,sizeof(uint64_t));
1246     stats->ins_cycles    = calloc(stats->nbases,sizeof(uint64_t));
1247     stats->del_cycles    = calloc(stats->nbases,sizeof(uint64_t));
1248     if ( targets )
1249         init_regions(stats, targets);
1250
1251     // Collect statistics
1252     if ( optind<argc )
1253     {
1254         // Collect stats in selected regions only
1255         bam_index_t *bam_idx = bam_index_load(bam_fname);
1256         if (bam_idx == 0)
1257             error("Random alignment retrieval only works for indexed BAM files.\n");
1258
1259         int i;
1260         for (i=optind; i<argc; i++) 
1261         {
1262             int tid, beg, end;
1263             bam_parse_region(stats->sam->header, argv[i], &tid, &beg, &end);
1264             if ( tid < 0 ) continue;
1265             bam_fetch(stats->sam->x.bam, bam_idx, tid, beg, end, stats, fetch_read);
1266         }
1267         bam_index_destroy(bam_idx);
1268     }
1269     else
1270     {
1271         // Stream through the entire BAM ignoring off-target regions if -t is given
1272         while (samread(sam,bam_line) >= 0) 
1273         {
1274             if ( stats->regions )
1275             {
1276                 if ( bam_line->core.tid >= stats->nregions || bam_line->core.tid<0 ) continue;
1277                 if ( !stats->is_sorted ) error("The BAM must be sorted in order for -t to work.\n");
1278
1279                 regions_t *reg = &stats->regions[bam_line->core.tid];
1280                 if ( reg->cpos==reg->npos ) continue;       // done for this chr
1281
1282                 // Find a matching interval or skip this read. No splicing of reads is done, no indels or soft clips considered, 
1283                 //  even small overlap is enough to include the read in the stats.
1284                 int i = reg->cpos;
1285                 while ( i<reg->npos && reg->pos[i].to<=bam_line->core.pos ) i++;
1286                 if ( i>=reg->npos ) { reg->cpos = reg->npos; continue; }
1287                 if ( bam_line->core.pos + bam_line->core.l_qseq + 1 < reg->pos[i].from ) continue;
1288                 reg->cpos = i;
1289             }
1290             collect_stats(bam_line,stats);
1291         }
1292     }
1293     round_buffer_flush(stats,-1);
1294
1295     output_stats(stats);
1296
1297     bam_destroy1(bam_line);
1298     samclose(stats->sam);
1299     if (stats->fai) fai_destroy(stats->fai);
1300     free(stats->cov_rbuf.buffer); free(stats->cov);
1301     free(stats->quals_1st); free(stats->quals_2nd); 
1302     free(stats->gc_1st); free(stats->gc_2nd);
1303     free(stats->isize_inward); free(stats->isize_outward); free(stats->isize_other);
1304     free(stats->gcd);
1305     free(stats->rseq_buf);
1306     free(stats->mpc_buf);
1307     free(stats->acgt_cycles);
1308     free(stats->read_lengths);
1309     free(stats->insertions);
1310     free(stats->deletions);
1311     free(stats->ins_cycles);
1312     free(stats->del_cycles);
1313     destroy_regions(stats);
1314     free(stats);
1315
1316         return 0;
1317 }
1318
1319
1320