2 Author: petr.danecek@sanger
3 gcc -Wall -Winline -g -O2 -I ~/git/samtools bamcheck.c -o bamcheck -lm -lz -L ~/git/samtools -lbam
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.
15 #define BAMCHECK_VERSION "2012-03-22"
17 #define _ISOC99_SOURCE
32 #define BWA_MIN_RDLEN 35
33 #define IS_PAIRED(bam) ((bam)->core.flag&BAM_FPAIRED && !((bam)->core.flag&BAM_FUNMAP) && !((bam)->core.flag&BAM_FMUNMAP))
34 #define IS_UNMAPPED(bam) ((bam)->core.flag&BAM_FUNMAP)
35 #define IS_REVERSE(bam) ((bam)->core.flag&BAM_FREVERSE)
36 #define IS_MATE_REVERSE(bam) ((bam)->core.flag&BAM_FMREVERSE)
37 #define IS_READ1(bam) ((bam)->core.flag&BAM_FREAD1)
38 #define IS_READ2(bam) ((bam)->core.flag&BAM_FREAD2)
39 #define IS_DUP(bam) ((bam)->core.flag&BAM_FDUP)
43 int32_t line_len, line_blen;
48 KHASH_MAP_INIT_STR(s, faidx1_t)
49 KHASH_MAP_INIT_STR(str, int)
64 // For coverage distribution, a simple pileup
73 typedef struct { uint32_t from, to; } pos_t;
84 int trim_qual; // bwa trim quality
85 int rmdup; // Exclude reads marked as duplicates from the stats
87 // Dimensions of the quality histogram holder (quals_1st,quals_2nd), GC content holder (gc_1st,gc_2nd),
88 // insert size histogram holder
89 int nquals; // The number of quality bins
90 int nbases; // The maximum sequence length the allocated array can hold
91 int nisize; // The maximum insert size that the allocated array can hold
92 int ngc; // The size of gc_1st and gc_2nd
93 int nindels; // The maximum indel length for indel distribution
95 // Arrays for the histogram data
96 uint64_t *quals_1st, *quals_2nd;
97 uint64_t *gc_1st, *gc_2nd;
98 uint64_t *isize_inward, *isize_outward, *isize_other;
99 uint64_t *acgt_cycles;
100 uint64_t *read_lengths;
101 uint64_t *insertions, *deletions;
102 uint64_t *ins_cycles, *del_cycles;
104 // The extremes encountered
105 int max_len; // Maximum read length
106 int max_qual; // Maximum quality
107 float isize_main_bulk; // There are always some unrealistically big insert sizes, report only the main part
112 uint64_t total_len_dup;
116 uint64_t nreads_unmapped;
117 uint64_t nreads_unpaired;
118 uint64_t nreads_paired;
120 uint64_t nbases_mapped;
121 uint64_t nbases_mapped_cigar;
122 uint64_t nbases_trimmed; // bwa trimmed bases
123 uint64_t nmismatches;
125 // GC-depth related data
126 uint32_t ngcd, igcd; // The maximum number of GC depth bins and index of the current bin
127 gc_depth_t *gcd; // The GC-depth bins holder
128 int gcd_bin_size; // The size of GC-depth bin
129 int32_t tid, gcd_pos; // Position of the current bin
130 int32_t pos; // Position of the last read
132 // Coverage distribution related data
133 int ncov; // The number of coverage bins
134 uint64_t *cov; // The coverage frequencies
135 int cov_min,cov_max,cov_step; // Minimum, maximum coverage and size of the coverage bins
136 round_buffer_t cov_rbuf; // Pileup round buffer
138 // Mismatches by read cycle
139 uint8_t *rseq_buf; // A buffer for reference sequence to check the mismatches against
140 int nref_seq; // The size of the buffer
141 int32_t rseq_pos; // The coordinate of the first base in the buffer
142 int32_t rseq_len; // The used part of the buffer
143 uint64_t *mpc_buf; // Mismatches per cycle
153 double sum_qual; // For calculating average quality value
154 samfile_t *sam; // Unused
155 faidx_t *fai; // Reference sequence for GC-depth graph
156 int argc; // Command line arguments to be printed on the output
161 void error(const char *format, ...);
163 // Coverage distribution methods
164 inline int coverage_idx(int min, int max, int n, int step, int depth)
172 return 1 + (depth - min) / step;
175 inline int round_buffer_lidx2ridx(int offset, int size, int64_t refpos, int64_t pos)
177 return (offset + (pos-refpos) % size) % size;
180 void round_buffer_flush(stats_t *stats, int64_t pos)
184 if ( pos==stats->cov_rbuf.pos )
187 int64_t new_pos = pos;
188 if ( pos==-1 || pos - stats->cov_rbuf.pos >= stats->cov_rbuf.size )
190 // Flush the whole buffer, but in sequential order,
191 pos = stats->cov_rbuf.pos + stats->cov_rbuf.size - 1;
194 if ( pos < stats->cov_rbuf.pos )
195 error("Expected coordinates in ascending order, got %ld after %ld\n", pos,stats->cov_rbuf.pos);
197 int ifrom = stats->cov_rbuf.start;
198 int ito = round_buffer_lidx2ridx(stats->cov_rbuf.start,stats->cov_rbuf.size,stats->cov_rbuf.pos,pos-1);
201 for (ibuf=ifrom; ibuf<stats->cov_rbuf.size; ibuf++)
203 if ( !stats->cov_rbuf.buffer[ibuf] )
205 idp = coverage_idx(stats->cov_min,stats->cov_max,stats->ncov,stats->cov_step,stats->cov_rbuf.buffer[ibuf]);
207 stats->cov_rbuf.buffer[ibuf] = 0;
211 for (ibuf=ifrom; ibuf<=ito; ibuf++)
213 if ( !stats->cov_rbuf.buffer[ibuf] )
215 idp = coverage_idx(stats->cov_min,stats->cov_max,stats->ncov,stats->cov_step,stats->cov_rbuf.buffer[ibuf]);
217 stats->cov_rbuf.buffer[ibuf] = 0;
219 stats->cov_rbuf.start = (new_pos==-1) ? 0 : round_buffer_lidx2ridx(stats->cov_rbuf.start,stats->cov_rbuf.size,stats->cov_rbuf.pos,pos);
220 stats->cov_rbuf.pos = new_pos;
223 void round_buffer_insert_read(round_buffer_t *rbuf, int64_t from, int64_t to)
225 if ( to-from >= rbuf->size )
226 error("The read length too big (%d), please increase the buffer length (currently %d)\n", to-from+1,rbuf->size);
227 if ( from < rbuf->pos )
228 error("The reads are not sorted (%ld comes after %ld).\n", from,rbuf->pos);
231 ifrom = round_buffer_lidx2ridx(rbuf->start,rbuf->size,rbuf->pos,from);
232 ito = round_buffer_lidx2ridx(rbuf->start,rbuf->size,rbuf->pos,to);
235 for (ibuf=ifrom; ibuf<rbuf->size; ibuf++)
236 rbuf->buffer[ibuf]++;
239 for (ibuf=ifrom; ibuf<=ito; ibuf++)
240 rbuf->buffer[ibuf]++;
243 // Calculate the number of bases in the read trimmed by BWA
244 int bwa_trim_read(int trim_qual, uint8_t *quals, int len, int reverse)
246 if ( len<BWA_MIN_RDLEN ) return 0;
248 // Although the name implies that the read cannot be trimmed to more than BWA_MIN_RDLEN,
249 // the calculation can in fact trim it to (BWA_MIN_RDLEN-1). (bwa_trim_read in bwa/bwaseqio.c).
250 int max_trimmed = len - BWA_MIN_RDLEN + 1;
251 int l, sum=0, max_sum=0, max_l=0;
253 for (l=0; l<max_trimmed; l++)
255 sum += trim_qual - quals[ reverse ? l : len-1-l ];
260 // This is the correct way, but bwa clips from some reason one base less
269 void count_indels(stats_t *stats,bam1_t *bam_line)
271 int is_fwd = IS_REVERSE(bam_line) ? 0 : 1;
274 int read_len = bam_line->core.l_qseq;
275 for (icig=0; icig<bam_line->core.n_cigar; icig++)
277 // Conversion from uint32_t to MIDNSHP
280 int cig = bam1_cigar(bam_line)[icig] & BAM_CIGAR_MASK;
281 int ncig = bam1_cigar(bam_line)[icig] >> BAM_CIGAR_SHIFT;
285 int idx = is_fwd ? icycle : read_len-icycle-1;
286 if ( idx >= stats->nbases ) error("FIXME: %d vs %d\n", idx,stats->nbases);
287 stats->ins_cycles[idx]++;
289 if ( ncig<=stats->nindels )
290 stats->insertions[ncig-1]++;
295 int idx = is_fwd ? icycle : read_len-icycle-1;
296 if ( idx >= stats->nbases ) error("FIXME: %d vs %d\n", idx,stats->nbases);
297 stats->del_cycles[idx]++;
298 if ( ncig<=stats->nindels )
299 stats->deletions[ncig-1]++;
306 void count_mismatches_per_cycle(stats_t *stats,bam1_t *bam_line)
308 int is_fwd = IS_REVERSE(bam_line) ? 0 : 1;
309 int icig,iread=0,icycle=0;
310 int iref = bam_line->core.pos - stats->rseq_pos;
311 int read_len = bam_line->core.l_qseq;
312 uint8_t *read = bam1_seq(bam_line);
313 uint8_t *quals = bam1_qual(bam_line);
314 uint64_t *mpc_buf = stats->mpc_buf;
315 for (icig=0; icig<bam_line->core.n_cigar; icig++)
317 // Conversion from uint32_t to MIDNSHP
320 int cig = bam1_cigar(bam_line)[icig] & BAM_CIGAR_MASK;
321 int ncig = bam1_cigar(bam_line)[icig] >> BAM_CIGAR_SHIFT;
336 // Soft-clips are present in the sequence, but the position of the read marks a start of non-clipped sequence
347 error("TODO: cigar %d, %s\n", cig,bam1_qname(bam_line));
349 if ( ncig+iref > stats->rseq_len )
350 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);
353 for (im=0; im<ncig; im++)
355 uint8_t cread = bam1_seqi(read,iread);
356 uint8_t cref = stats->rseq_buf[iref];
362 int idx = is_fwd ? icycle : read_len-icycle-1;
363 if ( idx>stats->max_len )
364 error("mpc: %d>%d\n",idx,stats->max_len);
365 idx = idx*stats->nquals;
366 if ( idx>=stats->nquals*stats->nbases )
367 error("FIXME: mpc_buf overflow\n");
370 else if ( cref && cread && cref!=cread )
372 uint8_t qual = quals[iread] + 1;
373 if ( qual>=stats->nquals )
374 error("TODO: quality too high %d>=%d\n", quals[iread],stats->nquals);
376 int idx = is_fwd ? icycle : read_len-icycle-1;
377 if ( idx>stats->max_len )
378 error("mpc: %d>%d\n",idx,stats->max_len);
380 idx = idx*stats->nquals + qual;
381 if ( idx>=stats->nquals*stats->nbases )
382 error("FIXME: mpc_buf overflow\n");
393 void read_ref_seq(stats_t *stats,int32_t tid,int32_t pos)
399 faidx_t *fai = stats->fai;
402 chr = stats->sam->header->target_name[tid];
404 // ID of the sequence name
405 iter = kh_get(s, h, chr);
406 if (iter == kh_end(h))
407 error("No such reference sequence [%s]?\n", chr);
408 val = kh_value(h, iter);
410 // Check the boundaries
412 error("Was the bam file mapped with the reference sequence supplied?"
413 " A read mapped beyond the end of the chromosome (%s:%d, chromosome length %d).\n", chr,pos,val.len);
414 int size = stats->nref_seq;
415 // The buffer extends beyond the chromosome end. Later the rest will be filled with N's.
416 if (size+pos > val.len) size = val.len-pos;
418 // Position the razf reader
419 razf_seek(fai->rz, val.offset + pos / val.line_blen * val.line_len + pos % val.line_blen, SEEK_SET);
421 uint8_t *ptr = stats->rseq_buf;
423 while ( nread<size && razf_read(fai->rz,&c,1) && !fai->rz->z_err )
428 // Conversion between uint8_t coding and ACGT
431 if ( c=='A' || c=='a' )
433 else if ( c=='C' || c=='c' )
435 else if ( c=='G' || c=='g' )
437 else if ( c=='T' || c=='t' )
444 if ( nread < stats->nref_seq )
446 memset(ptr,0, stats->nref_seq - nread);
447 nread = stats->nref_seq;
449 stats->rseq_len = nread;
450 stats->rseq_pos = pos;
454 float fai_gc_content(stats_t *stats)
457 int i,size = stats->rseq_len;
461 for (i=0; i<size; i++)
463 c = stats->rseq_buf[i];
469 else if ( c==1 || c==8 )
472 return count ? (float)gc/count : 0;
476 void realloc_buffers(stats_t *stats, int seq_len)
478 int n = 2*(1 + seq_len - stats->nbases) + stats->nbases;
480 stats->quals_1st = realloc(stats->quals_1st, n*stats->nquals*sizeof(uint64_t));
481 if ( !stats->quals_1st )
482 error("Could not realloc buffers, the sequence too long: %d (%ld)\n", seq_len,n*stats->nquals*sizeof(uint64_t));
483 memset(stats->quals_1st + stats->nbases*stats->nquals, 0, (n-stats->nbases)*stats->nquals*sizeof(uint64_t));
485 stats->quals_2nd = realloc(stats->quals_2nd, n*stats->nquals*sizeof(uint64_t));
486 if ( !stats->quals_2nd )
487 error("Could not realloc buffers, the sequence too long: %d (2x%ld)\n", seq_len,n*stats->nquals*sizeof(uint64_t));
488 memset(stats->quals_2nd + stats->nbases*stats->nquals, 0, (n-stats->nbases)*stats->nquals*sizeof(uint64_t));
490 if ( stats->mpc_buf )
492 stats->mpc_buf = realloc(stats->mpc_buf, n*stats->nquals*sizeof(uint64_t));
493 if ( !stats->mpc_buf )
494 error("Could not realloc buffers, the sequence too long: %d (%ld)\n", seq_len,n*stats->nquals*sizeof(uint64_t));
495 memset(stats->mpc_buf + stats->nbases*stats->nquals, 0, (n-stats->nbases)*stats->nquals*sizeof(uint64_t));
498 stats->acgt_cycles = realloc(stats->acgt_cycles, n*4*sizeof(uint64_t));
499 if ( !stats->acgt_cycles )
500 error("Could not realloc buffers, the sequence too long: %d (%ld)\n", seq_len,n*4*sizeof(uint64_t));
501 memset(stats->acgt_cycles + stats->nbases*4, 0, (n-stats->nbases)*4*sizeof(uint64_t));
503 stats->read_lengths = realloc(stats->read_lengths, n*sizeof(uint64_t));
504 if ( !stats->read_lengths )
505 error("Could not realloc buffers, the sequence too long: %d (%ld)\n", seq_len,n*sizeof(uint64_t));
506 memset(stats->read_lengths + stats->nbases, 0, (n-stats->nbases)*sizeof(uint64_t));
508 stats->insertions = realloc(stats->insertions, n*sizeof(uint64_t));
509 if ( !stats->insertions )
510 error("Could not realloc buffers, the sequence too long: %d (%ld)\n", seq_len,n*sizeof(uint64_t));
511 memset(stats->insertions + stats->nbases, 0, (n-stats->nbases)*sizeof(uint64_t));
513 stats->deletions = realloc(stats->deletions, n*sizeof(uint64_t));
514 if ( !stats->deletions )
515 error("Could not realloc buffers, the sequence too long: %d (%ld)\n", seq_len,n*sizeof(uint64_t));
516 memset(stats->deletions + stats->nbases, 0, (n-stats->nbases)*sizeof(uint64_t));
518 stats->ins_cycles = realloc(stats->ins_cycles, n*sizeof(uint64_t));
519 if ( !stats->ins_cycles )
520 error("Could not realloc buffers, the sequence too long: %d (%ld)\n", seq_len,n*sizeof(uint64_t));
521 memset(stats->ins_cycles + stats->nbases, 0, (n-stats->nbases)*sizeof(uint64_t));
523 stats->del_cycles = realloc(stats->del_cycles, n*sizeof(uint64_t));
524 if ( !stats->del_cycles )
525 error("Could not realloc buffers, the sequence too long: %d (%ld)\n", seq_len,n*sizeof(uint64_t));
526 memset(stats->del_cycles + stats->nbases, 0, (n-stats->nbases)*sizeof(uint64_t));
530 // Realloc the coverage distribution buffer
531 int *rbuffer = calloc(sizeof(int),seq_len*5);
532 n = stats->cov_rbuf.size-stats->cov_rbuf.start;
533 memcpy(rbuffer,stats->cov_rbuf.buffer+stats->cov_rbuf.start,n);
534 if ( stats->cov_rbuf.start>1 )
535 memcpy(rbuffer+n,stats->cov_rbuf.buffer,stats->cov_rbuf.start);
536 stats->cov_rbuf.start = 0;
537 free(stats->cov_rbuf.buffer);
538 stats->cov_rbuf.buffer = rbuffer;
539 stats->cov_rbuf.size = seq_len*5;
542 void collect_stats(bam1_t *bam_line, stats_t *stats)
544 if ( stats->rmdup && IS_DUP(bam_line) )
547 int seq_len = bam_line->core.l_qseq;
548 if ( !seq_len ) return;
549 if ( stats->filter_readlen!=-1 && seq_len!=stats->filter_readlen ) return;
550 if ( seq_len >= stats->nbases )
551 realloc_buffers(stats,seq_len);
552 if ( stats->max_len<seq_len )
553 stats->max_len = seq_len;
555 stats->read_lengths[seq_len]++;
557 // Count GC and ACGT per cycle
558 uint8_t base, *seq = bam1_seq(bam_line);
561 int reverse = IS_REVERSE(bam_line);
562 for (i=0; i<seq_len; i++)
564 // Conversion from uint8_t coding to ACGT
568 base = bam1_seqi(seq,i);
570 if ( base==1 || base==2 ) gc_count++;
571 else if ( base>2 ) base=3;
572 if ( 4*(reverse ? seq_len-i-1 : i) + base >= stats->nbases*4 )
573 error("FIXME: acgt_cycles\n");
574 stats->acgt_cycles[ 4*(reverse ? seq_len-i-1 : i) + base ]++;
576 int gc_idx_min = gc_count*(stats->ngc-1)/seq_len;
577 int gc_idx_max = (gc_count+1)*(stats->ngc-1)/seq_len;
578 if ( gc_idx_max >= stats->ngc ) gc_idx_max = stats->ngc - 1;
580 // Determine which array (1st or 2nd read) will these stats go to,
581 // trim low quality bases from end the same way BWA does,
584 uint8_t *bam_quals = bam1_qual(bam_line);
585 if ( bam_line->core.flag&BAM_FREAD2 )
587 quals = stats->quals_2nd;
589 for (i=gc_idx_min; i<gc_idx_max; i++)
594 quals = stats->quals_1st;
596 for (i=gc_idx_min; i<gc_idx_max; i++)
599 if ( stats->trim_qual>0 )
600 stats->nbases_trimmed += bwa_trim_read(stats->trim_qual, bam_quals, seq_len, reverse);
602 // Quality histogram and average quality
603 for (i=0; i<seq_len; i++)
605 uint8_t qual = bam_quals[ reverse ? seq_len-i-1 : i];
606 if ( qual>=stats->nquals )
607 error("TODO: quality too high %d>=%d\n", quals[i],stats->nquals);
608 if ( qual>stats->max_qual )
609 stats->max_qual = qual;
611 quals[ i*stats->nquals+qual ]++;
612 stats->sum_qual += qual;
615 // Look at the flags and increment appropriate counters (mapped, paired, etc)
616 if ( IS_UNMAPPED(bam_line) )
617 stats->nreads_unmapped++;
620 if ( !bam_line->core.qual )
623 count_indels(stats,bam_line);
625 // The insert size is tricky, because for long inserts the libraries are
626 // prepared differently and the pairs point in other direction. BWA does
627 // not set the paired flag for them. Similar thing is true also for 454
628 // reads. Therefore, do the insert size stats for all mapped reads.
629 int32_t isize = bam_line->core.isize;
630 if ( isize<0 ) isize = -isize;
631 if ( IS_PAIRED(bam_line) && isize!=0 )
633 stats->nreads_paired++;
634 if ( isize >= stats->nisize )
635 isize=stats->nisize-1;
637 int pos_fst = bam_line->core.mpos - bam_line->core.pos;
638 int is_fst = IS_READ1(bam_line) ? 1 : -1;
639 int is_fwd = IS_REVERSE(bam_line) ? -1 : 1;
640 int is_mfwd = IS_MATE_REVERSE(bam_line) ? -1 : 1;
642 if ( is_fwd*is_mfwd>0 )
643 stats->isize_other[isize]++;
644 else if ( is_fst*pos_fst>0 )
646 if ( is_fst*is_fwd>0 )
647 stats->isize_inward[isize]++;
649 stats->isize_outward[isize]++;
651 else if ( is_fst*pos_fst<0 )
653 if ( is_fst*is_fwd>0 )
654 stats->isize_outward[isize]++;
656 stats->isize_inward[isize]++;
660 stats->nreads_unpaired++;
662 // Number of mismatches
663 uint8_t *nm = bam_aux_get(bam_line,"NM");
665 stats->nmismatches += bam_aux2i(nm);
667 // Number of mapped bases from cigar
668 if ( bam_line->core.n_cigar == 0)
669 error("FIXME: mapped read with no cigar?\n");
670 int readlen = seq_len;
671 for (i=0; i<bam_line->core.n_cigar; i++)
673 // Conversion from uint32_t to MIDNSHP
676 if ( (bam1_cigar(bam_line)[i]&BAM_CIGAR_MASK)==0 || (bam1_cigar(bam_line)[i]&BAM_CIGAR_MASK)==1 )
677 stats->nbases_mapped_cigar += bam1_cigar(bam_line)[i]>>BAM_CIGAR_SHIFT;
679 if ( (bam1_cigar(bam_line)[i]&BAM_CIGAR_MASK)==2 )
680 readlen += bam1_cigar(bam_line)[i]>>BAM_CIGAR_SHIFT;
682 stats->nbases_mapped += seq_len;
684 if ( stats->tid==bam_line->core.tid && bam_line->core.pos<stats->pos )
685 stats->is_sorted = 0;
686 stats->pos = bam_line->core.pos;
688 if ( stats->is_sorted )
690 if ( stats->tid==-1 || stats->tid!=bam_line->core.tid )
691 round_buffer_flush(stats,-1);
693 // Mismatches per cycle and GC-depth graph
696 // First pass, new chromosome or sequence spanning beyond the end of the buffer
697 if ( stats->rseq_pos==-1 || stats->tid != bam_line->core.tid || stats->rseq_pos+stats->rseq_len < bam_line->core.pos+readlen )
699 read_ref_seq(stats,bam_line->core.tid,bam_line->core.pos);
701 // Get the reference GC content for this bin. Note that in the current implementation
702 // the GCD bins overlap by seq_len-1. Because the bin size is by default 20k and the
703 // expected read length only ~100bp, it shouldn't really matter.
704 stats->gcd_pos = bam_line->core.pos;
707 if ( stats->igcd >= stats->ngcd )
708 error("The genome too long or the GCD bin overlaps too big [%ud]\n", stats->igcd);
710 stats->gcd[ stats->igcd ].gc = fai_gc_content(stats);
712 count_mismatches_per_cycle(stats,bam_line);
714 else if ( stats->gcd_pos==-1 || stats->tid != bam_line->core.tid || bam_line->core.pos - stats->gcd_pos > stats->gcd_bin_size )
716 // First pass or a new chromosome
717 stats->tid = bam_line->core.tid;
718 stats->gcd_pos = bam_line->core.pos;
721 if ( stats->igcd >= stats->ngcd )
723 uint32_t n = 2*(1 + stats->ngcd);
724 stats->gcd = realloc(stats->gcd, n*sizeof(gc_depth_t));
726 error("Could not realloc GCD buffer, too many chromosomes or the genome too long?? [%u %u]\n", stats->ngcd,n);
727 memset(&(stats->gcd[stats->ngcd]),0,(n-stats->ngcd)*sizeof(gc_depth_t));
731 stats->gcd[ stats->igcd ].depth++;
732 // When no reference sequence is given, approximate the GC from the read (much shorter window, but otherwise OK)
734 stats->gcd[ stats->igcd ].gc += (float) gc_count / seq_len;
736 // Coverage distribution graph
737 round_buffer_flush(stats,bam_line->core.pos);
738 round_buffer_insert_read(&(stats->cov_rbuf),bam_line->core.pos,bam_line->core.pos+seq_len-1);
742 stats->total_len += seq_len;
743 if ( IS_DUP(bam_line) )
745 stats->total_len_dup += seq_len;
750 // Sort by GC and depth
751 #define GCD_t(x) ((gc_depth_t *)x)
752 static int gcd_cmp(const void *a, const void *b)
754 if ( GCD_t(a)->gc < GCD_t(b)->gc ) return -1;
755 if ( GCD_t(a)->gc > GCD_t(b)->gc ) return 1;
756 if ( GCD_t(a)->depth < GCD_t(b)->depth ) return -1;
757 if ( GCD_t(a)->depth > GCD_t(b)->depth ) return 1;
762 float gcd_percentile(gc_depth_t *gcd, int N, int p)
772 return gcd[N-1].depth;
775 return gcd[k-1].depth + d*(gcd[k].depth - gcd[k-1].depth);
778 void output_stats(stats_t *stats)
780 // Calculate average insert size and standard deviation (from the main bulk data only)
782 uint64_t nisize=0, nisize_inward=0, nisize_outward=0, nisize_other=0;
783 for (isize=1; isize<stats->nisize; isize++)
785 // Each pair was counted twice
786 stats->isize_inward[isize] *= 0.5;
787 stats->isize_outward[isize] *= 0.5;
788 stats->isize_other[isize] *= 0.5;
790 nisize_inward += stats->isize_inward[isize];
791 nisize_outward += stats->isize_outward[isize];
792 nisize_other += stats->isize_other[isize];
793 nisize += stats->isize_inward[isize] + stats->isize_outward[isize] + stats->isize_other[isize];
796 double bulk=0, avg_isize=0, sd_isize=0;
797 for (isize=1; isize<stats->nisize; isize++)
799 bulk += stats->isize_inward[isize] + stats->isize_outward[isize] + stats->isize_other[isize];
800 avg_isize += isize * (stats->isize_inward[isize] + stats->isize_outward[isize] + stats->isize_other[isize]);
802 if ( bulk/nisize > stats->isize_main_bulk )
809 avg_isize /= nisize ? nisize : 1;
810 for (isize=1; isize<ibulk; isize++)
811 sd_isize += (stats->isize_inward[isize] + stats->isize_outward[isize] + stats->isize_other[isize]) * (isize-avg_isize)*(isize-avg_isize) / nisize;
812 sd_isize = sqrt(sd_isize);
815 printf("# This file was produced by bamcheck (%s)\n",BAMCHECK_VERSION);
816 printf("# The command line was: %s",stats->argv[0]);
818 for (i=1; i<stats->argc; i++)
819 printf(" %s",stats->argv[i]);
821 printf("# Summary Numbers. Use `grep ^SN | cut -f 2-` to extract this part.\n");
822 printf("SN\tsequences:\t%ld\n", stats->nreads_1st+stats->nreads_2nd);
823 printf("SN\tis paired:\t%d\n", stats->nreads_1st&&stats->nreads_2nd ? 1 : 0);
824 printf("SN\tis sorted:\t%d\n", stats->is_sorted ? 1 : 0);
825 printf("SN\t1st fragments:\t%ld\n", stats->nreads_1st);
826 printf("SN\tlast fragments:\t%ld\n", stats->nreads_2nd);
827 printf("SN\treads mapped:\t%ld\n", stats->nreads_paired+stats->nreads_unpaired);
828 printf("SN\treads unmapped:\t%ld\n", stats->nreads_unmapped);
829 printf("SN\treads unpaired:\t%ld\n", stats->nreads_unpaired);
830 printf("SN\treads paired:\t%ld\n", stats->nreads_paired);
831 printf("SN\treads duplicated:\t%ld\n", stats->nreads_dup);
832 printf("SN\treads MQ0:\t%ld\n", stats->nreads_mq0);
833 printf("SN\ttotal length:\t%ld\n", stats->total_len);
834 printf("SN\tbases mapped:\t%ld\n", stats->nbases_mapped);
835 printf("SN\tbases mapped (cigar):\t%ld\n", stats->nbases_mapped_cigar);
836 printf("SN\tbases trimmed:\t%ld\n", stats->nbases_trimmed);
837 printf("SN\tbases duplicated:\t%ld\n", stats->total_len_dup);
838 printf("SN\tmismatches:\t%ld\n", stats->nmismatches);
839 printf("SN\terror rate:\t%e\n", (float)stats->nmismatches/stats->nbases_mapped_cigar);
840 float avg_read_length = (stats->nreads_1st+stats->nreads_2nd)?stats->total_len/(stats->nreads_1st+stats->nreads_2nd):0;
841 printf("SN\taverage length:\t%.0f\n", avg_read_length);
842 printf("SN\tmaximum length:\t%d\n", stats->max_len);
843 printf("SN\taverage quality:\t%.1f\n", stats->total_len?stats->sum_qual/stats->total_len:0);
844 printf("SN\tinsert size average:\t%.1f\n", avg_isize);
845 printf("SN\tinsert size standard deviation:\t%.1f\n", sd_isize);
846 printf("SN\tinward oriented pairs:\t%ld\n", nisize_inward);
847 printf("SN\toutward oriented pairs:\t%ld\n", nisize_outward);
848 printf("SN\tpairs with other orientation:\t%ld\n", nisize_other);
851 if ( stats->max_len<stats->nbases ) stats->max_len++;
852 if ( stats->max_qual+1<stats->nquals ) stats->max_qual++;
853 printf("# First Fragment Qualitites. Use `grep ^FFQ | cut -f 2-` to extract this part.\n");
854 printf("# Columns correspond to qualities and rows to cycles. First column is the cycle number.\n");
855 for (ibase=0; ibase<stats->max_len; ibase++)
857 printf("FFQ\t%d",ibase+1);
858 for (iqual=0; iqual<=stats->max_qual; iqual++)
860 printf("\t%ld", stats->quals_1st[ibase*stats->nquals+iqual]);
864 printf("# Last Fragment Qualitites. Use `grep ^LFQ | cut -f 2-` to extract this part.\n");
865 printf("# Columns correspond to qualities and rows to cycles. First column is the cycle number.\n");
866 for (ibase=0; ibase<stats->max_len; ibase++)
868 printf("LFQ\t%d",ibase+1);
869 for (iqual=0; iqual<=stats->max_qual; iqual++)
871 printf("\t%ld", stats->quals_2nd[ibase*stats->nquals+iqual]);
875 if ( stats->mpc_buf )
877 printf("# Mismatches per cycle and quality. Use `grep ^MPC | cut -f 2-` to extract this part.\n");
878 printf("# Columns correspond to qualities, rows to cycles. First column is the cycle number, second\n");
879 printf("# is the number of N's and the rest is the number of mismatches\n");
880 for (ibase=0; ibase<stats->max_len; ibase++)
882 printf("MPC\t%d",ibase+1);
883 for (iqual=0; iqual<=stats->max_qual; iqual++)
885 printf("\t%ld", stats->mpc_buf[ibase*stats->nquals+iqual]);
890 printf("# GC Content of first fragments. Use `grep ^GCF | cut -f 2-` to extract this part.\n");
892 for (ibase=0; ibase<stats->ngc; ibase++)
894 if ( stats->gc_1st[ibase]==stats->gc_1st[ibase_prev] ) continue;
895 printf("GCF\t%.2f\t%ld\n", (ibase+ibase_prev)*0.5*100./(stats->ngc-1),stats->gc_1st[ibase_prev]);
898 printf("# GC Content of last fragments. Use `grep ^GCL | cut -f 2-` to extract this part.\n");
900 for (ibase=0; ibase<stats->ngc; ibase++)
902 if ( stats->gc_2nd[ibase]==stats->gc_2nd[ibase_prev] ) continue;
903 printf("GCL\t%.2f\t%ld\n", (ibase+ibase_prev)*0.5*100./(stats->ngc-1),stats->gc_2nd[ibase_prev]);
906 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");
907 for (ibase=0; ibase<stats->max_len; ibase++)
909 uint64_t *ptr = &(stats->acgt_cycles[ibase*4]);
910 uint64_t sum = ptr[0]+ptr[1]+ptr[2]+ptr[3];
911 if ( ! sum ) continue;
912 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);
914 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");
915 for (isize=1; isize<ibulk; isize++)
916 printf("IS\t%d\t%ld\t%ld\t%ld\t%ld\n", isize,(stats->isize_inward[isize]+stats->isize_outward[isize]+stats->isize_other[isize]),
917 stats->isize_inward[isize],stats->isize_outward[isize],stats->isize_other[isize]);
919 printf("# Read lengths. Use `grep ^RL | cut -f 2-` to extract this part. The columns are: read length, count\n");
921 for (ilen=0; ilen<stats->max_len; ilen++)
923 if ( stats->read_lengths[ilen]>0 )
924 printf("RL\t%d\t%ld\n", ilen,stats->read_lengths[ilen]);
927 printf("# Indel distribution. Use `grep ^ID | cut -f 2-` to extract this part. The columns are: length, number of insertions, number of deletions\n");
928 for (ilen=0; ilen<stats->nindels; ilen++)
930 if ( stats->insertions[ilen]>0 || stats->deletions[ilen]>0 )
931 printf("ID\t%d\t%ld\t%ld\n", ilen+1,stats->insertions[ilen],stats->deletions[ilen]);
934 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");
935 for (ilen=0; ilen<stats->nbases; ilen++)
937 if ( stats->ins_cycles[ilen]>0 || stats->del_cycles[ilen]>0 )
938 printf("IC\t%d\t%ld\t%ld\n", ilen+1,stats->ins_cycles[ilen],stats->del_cycles[ilen]);
941 printf("# Coverage distribution. Use `grep ^COV | cut -f 2-` to extract this part.\n");
942 printf("COV\t[<%d]\t%d\t%ld\n",stats->cov_min,stats->cov_min-1,stats->cov[0]);
944 for (icov=1; icov<stats->ncov-1; icov++)
945 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,stats->cov[icov]);
946 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,stats->cov[stats->ncov-1]);
949 // Calculate average GC content, then sort by GC and depth
950 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");
952 for (igcd=0; igcd<stats->igcd; igcd++)
955 stats->gcd[igcd].gc = round(100. * stats->gcd[igcd].gc);
957 if ( stats->gcd[igcd].depth )
958 stats->gcd[igcd].gc = round(100. * stats->gcd[igcd].gc / stats->gcd[igcd].depth);
960 qsort(stats->gcd, stats->igcd+1, sizeof(gc_depth_t), gcd_cmp);
962 while ( igcd < stats->igcd )
964 // Calculate percentiles (10,25,50,75,90th) for the current GC content and print
965 uint32_t nbins=0, itmp=igcd;
966 float gc = stats->gcd[igcd].gc;
967 while ( itmp<stats->igcd && fabs(stats->gcd[itmp].gc-gc)<0.1 )
972 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),
973 gcd_percentile(&(stats->gcd[igcd]),nbins,10) *avg_read_length/stats->gcd_bin_size,
974 gcd_percentile(&(stats->gcd[igcd]),nbins,25) *avg_read_length/stats->gcd_bin_size,
975 gcd_percentile(&(stats->gcd[igcd]),nbins,50) *avg_read_length/stats->gcd_bin_size,
976 gcd_percentile(&(stats->gcd[igcd]),nbins,75) *avg_read_length/stats->gcd_bin_size,
977 gcd_percentile(&(stats->gcd[igcd]),nbins,90) *avg_read_length/stats->gcd_bin_size
983 void bam_init_header_hash(bam_header_t *header);
985 void init_regions(stats_t *stats, char *file)
988 khash_t(str) *header_hash;
990 bam_init_header_hash(stats->sam->header);
991 header_hash = (khash_t(str)*)stats->sam->header->hash;
993 FILE *fp = fopen(file,"r");
994 if ( !fp ) error("%s: %s\n",file,strerror(errno));
1000 int prev_tid=-1, prev_pos=-1;
1001 while ((nread = getline(&line, &len, fp)) != -1)
1003 if ( line[0] == '#' ) continue;
1006 while ( i<nread && !isspace(line[i]) ) i++;
1007 if ( i>=nread ) error("Could not parse the file: %s\n", file);
1010 iter = kh_get(str, header_hash, line);
1011 int tid = kh_val(header_hash, iter);
1012 if ( iter == kh_end(header_hash) )
1015 fprintf(stderr,"Warning: Some sequences not present in the BAM (%s)\n", line);
1020 if ( tid >= stats->nregions )
1022 stats->regions = realloc(stats->regions,sizeof(regions_t)*(stats->nregions+100));
1024 for (j=stats->nregions; j<stats->nregions+100; j++)
1026 stats->regions[j].npos = stats->regions[j].mpos = stats->regions[j].cpos = 0;
1027 stats->regions[j].pos = NULL;
1029 stats->nregions += 100;
1031 int npos = stats->regions[tid].npos;
1032 if ( npos >= stats->regions[tid].mpos )
1034 stats->regions[tid].mpos += 1000;
1035 stats->regions[tid].pos = realloc(stats->regions[tid].pos,sizeof(pos_t)*stats->regions[tid].mpos);
1038 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");
1039 if ( prev_tid==-1 || prev_tid!=tid )
1042 prev_pos = stats->regions[tid].pos[npos].from;
1044 if ( prev_pos>stats->regions[tid].pos[npos].from )
1045 error("The positions are not in chromosomal order (%s:%d comes after %d)\n", line,stats->regions[tid].pos[npos].from,prev_pos);
1046 stats->regions[tid].npos++;
1048 if (line) free(line);
1052 void destroy_regions(stats_t *stats)
1055 for (i=0; i<stats->nregions; i++)
1057 if ( !stats->regions[i].mpos ) continue;
1058 free(stats->regions[i].pos);
1060 if ( stats->regions ) free(stats->regions);
1063 static int fetch_read(const bam1_t *bam_line, void *data)
1065 collect_stats((bam1_t*)bam_line,(stats_t*)data);
1070 void error(const char *format, ...)
1074 printf("Version: %s\n", BAMCHECK_VERSION);
1075 printf("About: The program collects statistics from BAM files. The output can be visualized using plot-bamcheck.\n");
1076 printf("Usage: bamcheck [OPTIONS] file.bam\n");
1077 printf(" bamcheck [OPTIONS] file.bam chr:from-to\n");
1078 printf("Options:\n");
1079 printf(" -c, --coverage <int>,<int>,<int> Coverage distribution min,max,step [1,1000,1]\n");
1080 printf(" -d, --remove-dups Exlude from statistics reads marked as duplicates\n");
1081 printf(" -h, --help This help message\n");
1082 printf(" -i, --insert-size <int> Maximum insert size [8000]\n");
1083 printf(" -l, --read-length <int> Include in the statistics only reads with the given read length []\n");
1084 printf(" -m, --most-inserts <float> Report only the main part of inserts [0.99]\n");
1085 printf(" -q, --trim-quality <int> The BWA trimming parameter [0]\n");
1086 printf(" -r, --ref-seq <file> Reference sequence (required for GC-depth calculation).\n");
1087 printf(" -t, --target-regions <file> Do stats in these regions only. Tab-delimited file chr,from,to, 1-based, inclusive.\n");
1088 printf(" -s, --sam Input is SAM\n");
1094 va_start(ap, format);
1095 vfprintf(stderr, format, ap);
1101 int main(int argc, char *argv[])
1103 char *targets = NULL;
1104 char *bam_fname = NULL;
1105 samfile_t *sam = NULL;
1108 stats_t *stats = calloc(1,sizeof(stats_t));
1111 stats->nbases = 300;
1112 stats->nisize = 8000;
1113 stats->max_len = 30;
1114 stats->max_qual = 40;
1115 stats->isize_main_bulk = 0.99; // There are always outliers at the far end
1116 stats->gcd_bin_size = 20000;
1117 stats->ngcd = 3e5; // 300k of 20k bins is enough to hold a genome 6Gbp big
1118 stats->nref_seq = stats->gcd_bin_size;
1119 stats->rseq_pos = -1;
1120 stats->tid = stats->gcd_pos = -1;
1121 stats->is_sorted = 1;
1123 stats->cov_max = 1000;
1124 stats->cov_step = 1;
1127 stats->filter_readlen = -1;
1128 stats->nindels = stats->nbases;
1130 strcpy(in_mode, "rb");
1132 static struct option loptions[] =
1135 {"remove-dups",0,0,'d'},
1137 {"ref-seq",0,0,'r'},
1138 {"coverage",0,0,'c'},
1139 {"read-length",0,0,'l'},
1140 {"insert-size",0,0,'i'},
1141 {"most-inserts",0,0,'m'},
1142 {"trim-quality",0,0,'q'},
1143 {"target-regions",0,0,'t'},
1147 while ( (opt=getopt_long(argc,argv,"?hdsr:c:l:i:t:m:q:",loptions,NULL))>0 )
1151 case 'd': stats->rmdup=1; break;
1152 case 's': strcpy(in_mode, "r"); break;
1153 case 'r': stats->fai = fai_load(optarg);
1155 error("Could not load faidx: %s\n", optarg);
1157 case 'c': if ( sscanf(optarg,"%d,%d,%d",&stats->cov_min,&stats->cov_max,&stats->cov_step)!= 3 )
1158 error("Unable to parse -c %s\n", optarg);
1160 case 'l': stats->filter_readlen = atoi(optarg); break;
1161 case 'i': stats->nisize = atoi(optarg); break;
1162 case 'm': stats->isize_main_bulk = atof(optarg); break;
1163 case 'q': stats->trim_qual = atoi(optarg); break;
1164 case 't': targets = optarg; break;
1166 case 'h': error(NULL);
1167 default: error("Unknown argument: %s\n", optarg);
1171 bam_fname = argv[optind++];
1175 if ( isatty(fileno((FILE *)stdin)) )
1181 // .. coverage bins and round buffer
1182 if ( stats->cov_step > stats->cov_max - stats->cov_min + 1 )
1184 stats->cov_step = stats->cov_max - stats->cov_min;
1185 if ( stats->cov_step <= 0 )
1186 stats->cov_step = 1;
1188 stats->ncov = 3 + (stats->cov_max-stats->cov_min) / stats->cov_step;
1189 stats->cov_max = stats->cov_min + ((stats->cov_max-stats->cov_min)/stats->cov_step +1)*stats->cov_step - 1;
1190 stats->cov = calloc(sizeof(uint64_t),stats->ncov);
1191 stats->cov_rbuf.size = stats->nbases*5;
1192 stats->cov_rbuf.buffer = calloc(sizeof(int32_t),stats->cov_rbuf.size);
1194 if ((sam = samopen(bam_fname, in_mode, NULL)) == 0)
1195 error("Failed to open: %s\n", bam_fname);
1197 bam1_t *bam_line = bam_init1();
1199 stats->quals_1st = calloc(stats->nquals*stats->nbases,sizeof(uint64_t));
1200 stats->quals_2nd = calloc(stats->nquals*stats->nbases,sizeof(uint64_t));
1201 stats->gc_1st = calloc(stats->ngc,sizeof(uint64_t));
1202 stats->gc_2nd = calloc(stats->ngc,sizeof(uint64_t));
1203 stats->isize_inward = calloc(stats->nisize,sizeof(uint64_t));
1204 stats->isize_outward = calloc(stats->nisize,sizeof(uint64_t));
1205 stats->isize_other = calloc(stats->nisize,sizeof(uint64_t));
1206 stats->gcd = calloc(stats->ngcd,sizeof(gc_depth_t));
1207 stats->rseq_buf = calloc(stats->nref_seq,sizeof(uint8_t));
1208 stats->mpc_buf = stats->fai ? calloc(stats->nquals*stats->nbases,sizeof(uint64_t)) : NULL;
1209 stats->acgt_cycles = calloc(4*stats->nbases,sizeof(uint64_t));
1210 stats->read_lengths = calloc(stats->nbases,sizeof(uint64_t));
1211 stats->insertions = calloc(stats->nbases,sizeof(uint64_t));
1212 stats->deletions = calloc(stats->nbases,sizeof(uint64_t));
1213 stats->ins_cycles = calloc(stats->nbases,sizeof(uint64_t));
1214 stats->del_cycles = calloc(stats->nbases,sizeof(uint64_t));
1216 init_regions(stats, targets);
1218 // Collect statistics
1221 // Collect stats in selected regions only
1222 bam_index_t *bam_idx = bam_index_load(bam_fname);
1224 error("Random alignment retrieval only works for indexed BAM files.\n");
1227 for (i=optind; i<argc; i++)
1230 bam_parse_region(stats->sam->header, argv[i], &tid, &beg, &end);
1231 if ( tid < 0 ) continue;
1232 bam_fetch(stats->sam->x.bam, bam_idx, tid, beg, end, stats, fetch_read);
1234 bam_index_destroy(bam_idx);
1238 // Stream through the entire BAM ignoring off-target regions if -t is given
1239 while (samread(sam,bam_line) >= 0)
1241 if ( stats->regions )
1243 if ( bam_line->core.tid >= stats->nregions ) continue;
1244 if ( !stats->is_sorted ) error("The BAM must be sorted in order for -t to work.\n");
1246 regions_t *reg = &stats->regions[bam_line->core.tid];
1247 if ( reg->cpos==reg->npos ) continue; // done for this chr
1249 // Find a matching interval or skip this read. No splicing of reads is done, no indels or soft clips considered,
1250 // even small overlap is enough to include the read in the stats.
1252 while ( i<reg->npos && reg->pos[i].to<=bam_line->core.pos ) i++;
1253 if ( i>=reg->npos ) { reg->cpos = reg->npos; continue; }
1254 if ( bam_line->core.pos + bam_line->core.l_qseq + 1 < reg->pos[i].from ) continue;
1257 collect_stats(bam_line,stats);
1260 round_buffer_flush(stats,-1);
1262 output_stats(stats);
1264 bam_destroy1(bam_line);
1265 samclose(stats->sam);
1266 if (stats->fai) fai_destroy(stats->fai);
1267 free(stats->cov_rbuf.buffer); free(stats->cov);
1268 free(stats->quals_1st); free(stats->quals_2nd);
1269 free(stats->gc_1st); free(stats->gc_2nd);
1270 free(stats->isize_inward); free(stats->isize_outward); free(stats->isize_other);
1272 free(stats->rseq_buf);
1273 free(stats->mpc_buf);
1274 free(stats->acgt_cycles);
1275 free(stats->read_lengths);
1276 free(stats->insertions);
1277 free(stats->deletions);
1278 free(stats->ins_cycles);
1279 free(stats->del_cycles);
1280 destroy_regions(stats);