Merge commit 'upstream/0.1.17'
[samtools.git] / sam_view.c
1 #include <stdlib.h>
2 #include <string.h>
3 #include <stdio.h>
4 #include <unistd.h>
5 #include <math.h>
6 #include "sam_header.h"
7 #include "sam.h"
8 #include "faidx.h"
9 #include "kstring.h"
10 #include "khash.h"
11 KHASH_SET_INIT_STR(rg)
12
13 // When counting records instead of printing them,
14 // data passed to the bam_fetch callback is encapsulated in this struct.
15 typedef struct {
16         bam_header_t *header;
17         int *count;
18 } count_func_data_t;
19
20 typedef khash_t(rg) *rghash_t;
21
22 static rghash_t g_rghash = 0;
23 static int g_min_mapQ = 0, g_flag_on = 0, g_flag_off = 0;
24 static char *g_library, *g_rg;
25 static void *g_bed;
26
27 void *bed_read(const char *fn);
28 void bed_destroy(void *_h);
29 int bed_overlap(const void *_h, const char *chr, int beg, int end);
30
31 static inline int __g_skip_aln(const bam_header_t *h, const bam1_t *b)
32 {
33         if (b->core.qual < g_min_mapQ || ((b->core.flag & g_flag_on) != g_flag_on) || (b->core.flag & g_flag_off))
34                 return 1;
35         if (g_bed && b->core.tid >= 0 && !bed_overlap(g_bed, h->target_name[b->core.tid], b->core.pos, bam_calend(&b->core, bam1_cigar(b))))
36                 return 1;
37         if (g_rg || g_rghash) {
38                 uint8_t *s = bam_aux_get(b, "RG");
39                 if (s) {
40                         if (g_rg) return (strcmp(g_rg, (char*)(s + 1)) == 0)? 0 : 1;
41                         if (g_rghash) {
42                                 khint_t k = kh_get(rg, g_rghash, (char*)(s + 1));
43                                 return (k != kh_end(g_rghash))? 0 : 1;
44                         }
45                 }
46         }
47         if (g_library) {
48                 const char *p = bam_get_library((bam_header_t*)h, b);
49                 return (p && strcmp(p, g_library) == 0)? 0 : 1;
50         }
51         return 0;
52 }
53
54 static char *drop_rg(char *hdtxt, rghash_t h, int *len)
55 {
56         char *p = hdtxt, *q, *r, *s;
57         kstring_t str;
58         memset(&str, 0, sizeof(kstring_t));
59         while (1) {
60                 int toprint = 0;
61                 q = strchr(p, '\n');
62                 if (q == 0) q = p + strlen(p);
63                 if (q - p < 3) break; // the line is too short; then stop
64                 if (strncmp(p, "@RG\t", 4) == 0) {
65                         int c;
66                         khint_t k;
67                         if ((r = strstr(p, "\tID:")) != 0) {
68                                 r += 4;
69                                 for (s = r; *s != '\0' && *s != '\n' && *s != '\t'; ++s);
70                                 c = *s; *s = '\0';
71                                 k = kh_get(rg, h, r);
72                                 *s = c;
73                                 if (k != kh_end(h)) toprint = 1;
74                         }
75                 } else toprint = 1;
76                 if (toprint) {
77                         kputsn(p, q - p, &str); kputc('\n', &str);
78                 }
79                 p = q + 1;
80         }
81         *len = str.l;
82         return str.s;
83 }
84
85 // callback function for bam_fetch() that prints nonskipped records
86 static int view_func(const bam1_t *b, void *data)
87 {
88         if (!__g_skip_aln(((samfile_t*)data)->header, b))
89                 samwrite((samfile_t*)data, b);
90         return 0;
91 }
92
93 // callback function for bam_fetch() that counts nonskipped records
94 static int count_func(const bam1_t *b, void *data)
95 {
96         if (!__g_skip_aln(((count_func_data_t*)data)->header, b)) {
97                 (*((count_func_data_t*)data)->count)++;
98         }
99         return 0;
100 }
101
102 static int usage(int is_long_help);
103
104 int main_samview(int argc, char *argv[])
105 {
106         int c, is_header = 0, is_header_only = 0, is_bamin = 1, ret = 0, compress_level = -1, is_bamout = 0, is_count = 0;
107         int of_type = BAM_OFDEC, is_long_help = 0;
108         int count = 0;
109         samfile_t *in = 0, *out = 0;
110         char in_mode[5], out_mode[5], *fn_out = 0, *fn_list = 0, *fn_ref = 0, *fn_rg = 0;
111
112         /* parse command-line options */
113         strcpy(in_mode, "r"); strcpy(out_mode, "w");
114         while ((c = getopt(argc, argv, "Sbct:h1Ho:q:f:F:ul:r:xX?T:R:L:")) >= 0) {
115                 switch (c) {
116                 case 'c': is_count = 1; break;
117                 case 'S': is_bamin = 0; break;
118                 case 'b': is_bamout = 1; break;
119                 case 't': fn_list = strdup(optarg); is_bamin = 0; break;
120                 case 'h': is_header = 1; break;
121                 case 'H': is_header_only = 1; break;
122                 case 'o': fn_out = strdup(optarg); break;
123                 case 'f': g_flag_on = strtol(optarg, 0, 0); break;
124                 case 'F': g_flag_off = strtol(optarg, 0, 0); break;
125                 case 'q': g_min_mapQ = atoi(optarg); break;
126                 case 'u': compress_level = 0; break;
127                 case '1': compress_level = 1; break;
128                 case 'l': g_library = strdup(optarg); break;
129                 case 'L': g_bed = bed_read(optarg); break;
130                 case 'r': g_rg = strdup(optarg); break;
131                 case 'R': fn_rg = strdup(optarg); break;
132                 case 'x': of_type = BAM_OFHEX; break;
133                 case 'X': of_type = BAM_OFSTR; break;
134                 case '?': is_long_help = 1; break;
135                 case 'T': fn_ref = strdup(optarg); is_bamin = 0; break;
136                 default: return usage(is_long_help);
137                 }
138         }
139         if (compress_level >= 0) is_bamout = 1;
140         if (is_header_only) is_header = 1;
141         if (is_bamout) strcat(out_mode, "b");
142         else {
143                 if (of_type == BAM_OFHEX) strcat(out_mode, "x");
144                 else if (of_type == BAM_OFSTR) strcat(out_mode, "X");
145         }
146         if (is_bamin) strcat(in_mode, "b");
147         if (is_header) strcat(out_mode, "h");
148         if (compress_level >= 0) {
149                 char tmp[2];
150                 tmp[0] = compress_level + '0'; tmp[1] = '\0';
151                 strcat(out_mode, tmp);
152         }
153         if (argc == optind) return usage(is_long_help); // potential memory leak...
154
155         // read the list of read groups
156         if (fn_rg) {
157                 FILE *fp_rg;
158                 char buf[1024];
159                 int ret;
160                 g_rghash = kh_init(rg);
161                 fp_rg = fopen(fn_rg, "r");
162                 while (!feof(fp_rg) && fscanf(fp_rg, "%s", buf) > 0) // this is not a good style, but bear me...
163                         kh_put(rg, g_rghash, strdup(buf), &ret); // we'd better check duplicates...
164                 fclose(fp_rg);
165         }
166
167         // generate the fn_list if necessary
168         if (fn_list == 0 && fn_ref) fn_list = samfaipath(fn_ref);
169         // open file handlers
170         if ((in = samopen(argv[optind], in_mode, fn_list)) == 0) {
171                 fprintf(stderr, "[main_samview] fail to open \"%s\" for reading.\n", argv[optind]);
172                 ret = 1;
173                 goto view_end;
174         }
175         if (in->header == 0) {
176                 fprintf(stderr, "[main_samview] fail to read the header from \"%s\".\n", argv[optind]);
177                 ret = 1;
178                 goto view_end;
179         }
180         if (g_rghash) { // FIXME: I do not know what "bam_header_t::n_text" is for...
181                 char *tmp;
182                 int l;
183                 tmp = drop_rg(in->header->text, g_rghash, &l);
184                 free(in->header->text);
185                 in->header->text = tmp;
186                 in->header->l_text = l;
187         }
188         if (!is_count && (out = samopen(fn_out? fn_out : "-", out_mode, in->header)) == 0) {
189                 fprintf(stderr, "[main_samview] fail to open \"%s\" for writing.\n", fn_out? fn_out : "standard output");
190                 ret = 1;
191                 goto view_end;
192         }
193         if (is_header_only) goto view_end; // no need to print alignments
194
195         if (argc == optind + 1) { // convert/print the entire file
196                 bam1_t *b = bam_init1();
197                 int r;
198                 while ((r = samread(in, b)) >= 0) { // read one alignment from `in'
199                         if (!__g_skip_aln(in->header, b)) {
200                                 if (!is_count) samwrite(out, b); // write the alignment to `out'
201                                 count++;
202                         }
203                 }
204                 if (r < -1) {
205                         fprintf(stderr, "[main_samview] truncated file.\n");
206                         ret = 1;
207                 }
208                 bam_destroy1(b);
209         } else { // retrieve alignments in specified regions
210                 int i;
211                 bam_index_t *idx = 0;
212                 if (is_bamin) idx = bam_index_load(argv[optind]); // load BAM index
213                 if (idx == 0) { // index is unavailable
214                         fprintf(stderr, "[main_samview] random alignment retrieval only works for indexed BAM files.\n");
215                         ret = 1;
216                         goto view_end;
217                 }
218                 for (i = optind + 1; i < argc; ++i) {
219                         int tid, beg, end, result;
220                         bam_parse_region(in->header, argv[i], &tid, &beg, &end); // parse a region in the format like `chr2:100-200'
221                         if (tid < 0) { // reference name is not found
222                                 fprintf(stderr, "[main_samview] region \"%s\" specifies an unknown reference name. Continue anyway.\n", argv[i]);
223                                 continue;
224                         }
225                         // fetch alignments
226                         if (is_count) {
227                                 count_func_data_t count_data = { in->header, &count };
228                                 result = bam_fetch(in->x.bam, idx, tid, beg, end, &count_data, count_func);
229                         } else
230                                 result = bam_fetch(in->x.bam, idx, tid, beg, end, out, view_func);
231                         if (result < 0) {
232                                 fprintf(stderr, "[main_samview] retrieval of region \"%s\" failed due to truncated file or corrupt BAM index file\n", argv[i]);
233                                 ret = 1;
234                                 break;
235                         }
236                 }
237                 bam_index_destroy(idx); // destroy the BAM index
238         }
239
240 view_end:
241         if (is_count && ret == 0) {
242                 printf("%d\n", count);
243         }
244         // close files, free and return
245         free(fn_list); free(fn_ref); free(fn_out); free(g_library); free(g_rg); free(fn_rg);
246         if (g_bed) bed_destroy(g_bed);
247         if (g_rghash) {
248                 khint_t k;
249                 for (k = 0; k < kh_end(g_rghash); ++k)
250                         if (kh_exist(g_rghash, k)) free((char*)kh_key(g_rghash, k));
251                 kh_destroy(rg, g_rghash);
252         }
253         samclose(in);
254         if (!is_count)
255                 samclose(out);
256         return ret;
257 }
258
259 static int usage(int is_long_help)
260 {
261         fprintf(stderr, "\n");
262         fprintf(stderr, "Usage:   samtools view [options] <in.bam>|<in.sam> [region1 [...]]\n\n");
263         fprintf(stderr, "Options: -b       output BAM\n");
264         fprintf(stderr, "         -h       print header for the SAM output\n");
265         fprintf(stderr, "         -H       print header only (no alignments)\n");
266         fprintf(stderr, "         -S       input is SAM\n");
267         fprintf(stderr, "         -u       uncompressed BAM output (force -b)\n");
268         fprintf(stderr, "         -1       fast compression (force -b)\n");
269         fprintf(stderr, "         -x       output FLAG in HEX (samtools-C specific)\n");
270         fprintf(stderr, "         -X       output FLAG in string (samtools-C specific)\n");
271         fprintf(stderr, "         -c       print only the count of matching records\n");
272         fprintf(stderr, "         -L FILE  output alignments overlapping the input BED FILE [null]\n");
273         fprintf(stderr, "         -t FILE  list of reference names and lengths (force -S) [null]\n");
274         fprintf(stderr, "         -T FILE  reference sequence file (force -S) [null]\n");
275         fprintf(stderr, "         -o FILE  output file name [stdout]\n");
276         fprintf(stderr, "         -R FILE  list of read groups to be outputted [null]\n");
277         fprintf(stderr, "         -f INT   required flag, 0 for unset [0]\n");
278         fprintf(stderr, "         -F INT   filtering flag, 0 for unset [0]\n");
279         fprintf(stderr, "         -q INT   minimum mapping quality [0]\n");
280         fprintf(stderr, "         -l STR   only output reads in library STR [null]\n");
281         fprintf(stderr, "         -r STR   only output reads in read group STR [null]\n");
282         fprintf(stderr, "         -?       longer help\n");
283         fprintf(stderr, "\n");
284         if (is_long_help)
285                 fprintf(stderr, "Notes:\n\
286 \n\
287   1. By default, this command assumes the file on the command line is in\n\
288      the BAM format and it prints the alignments in SAM. If `-t' is\n\
289      applied, the input file is assumed to be in the SAM format. The\n\
290      file supplied with `-t' is SPACE/TAB delimited with the first two\n\
291      fields of each line consisting of the reference name and the\n\
292      corresponding sequence length. The `.fai' file generated by `faidx'\n\
293      can be used here. This file may be empty if reads are unaligned.\n\
294 \n\
295   2. SAM->BAM conversion: `samtools view -bT ref.fa in.sam.gz'.\n\
296 \n\
297   3. BAM->SAM conversion: `samtools view in.bam'.\n\
298 \n\
299   4. A region should be presented in one of the following formats:\n\
300      `chr1', `chr2:1,000' and `chr3:1000-2,000'. When a region is\n\
301      specified, the input alignment file must be an indexed BAM file.\n\
302 \n\
303   5. Option `-u' is preferred over `-b' when the output is piped to\n\
304      another samtools command.\n\
305 \n\
306   6. In a string FLAG, each character represents one bit with\n\
307      p=0x1 (paired), P=0x2 (properly paired), u=0x4 (unmapped),\n\
308      U=0x8 (mate unmapped), r=0x10 (reverse), R=0x20 (mate reverse)\n\
309      1=0x40 (first), 2=0x80 (second), s=0x100 (not primary), \n\
310      f=0x200 (failure) and d=0x400 (duplicate). Note that `-x' and\n\
311      `-X' are samtools-C specific. Picard and older samtools do not\n\
312      support HEX or string flags.\n\
313 \n");
314         return 1;
315 }
316
317 int main_import(int argc, char *argv[])
318 {
319         int argc2, ret;
320         char **argv2;
321         if (argc != 4) {
322                 fprintf(stderr, "Usage: bamtk import <in.ref_list> <in.sam> <out.bam>\n");
323                 return 1;
324         }
325         argc2 = 6;
326         argv2 = calloc(6, sizeof(char*));
327         argv2[0] = "import", argv2[1] = "-o", argv2[2] = argv[3], argv2[3] = "-bt", argv2[4] = argv[1], argv2[5] = argv[2];
328         ret = main_samview(argc2, argv2);
329         free(argv2);
330         return ret;
331 }
332
333 int8_t seq_comp_table[16] = { 0, 8, 4, 12, 2, 10, 9, 14, 1, 6, 5, 13, 3, 11, 7, 15 };
334
335 int main_bam2fq(int argc, char *argv[])
336 {
337         bamFile fp;
338         bam_header_t *h;
339         bam1_t *b;
340         int8_t *buf;
341         int max_buf;
342         if (argc == 1) {
343                 fprintf(stderr, "Usage: samtools bam2fq <in.bam>\n");
344                 return 1;
345         }
346         fp = strcmp(argv[1], "-")? bam_open(argv[1], "r") : bam_dopen(fileno(stdin), "r");
347         if (fp == 0) return 1;
348         h = bam_header_read(fp);
349         b = bam_init1();
350         buf = 0;
351         max_buf = 0;
352         while (bam_read1(fp, b) >= 0) {
353                 int i, qlen = b->core.l_qseq;
354                 uint8_t *seq;
355                 putchar('@'); fputs(bam1_qname(b), stdout);
356                 if ((b->core.flag & 0x40) && !(b->core.flag & 0x80)) puts("/1");
357                 else if ((b->core.flag & 0x80) && !(b->core.flag & 0x40)) puts("/2");
358                 else putchar('\n');
359                 if (max_buf < qlen + 1) {
360                         max_buf = qlen + 1;
361                         kroundup32(max_buf);
362                         buf = realloc(buf, max_buf);
363                 }
364                 buf[qlen] = 0;
365                 seq = bam1_seq(b);
366                 for (i = 0; i < qlen; ++i)
367                         buf[i] = bam1_seqi(seq, i);
368                 if (b->core.flag & 16) { // reverse complement
369                         for (i = 0; i < qlen>>1; ++i) {
370                                 int8_t t = seq_comp_table[buf[qlen - 1 - i]];
371                                 buf[qlen - 1 - i] = seq_comp_table[buf[i]];
372                                 buf[i] = t;
373                         }
374                         if (qlen&1) buf[i] = seq_comp_table[buf[i]];
375                 }
376                 for (i = 0; i < qlen; ++i)
377                         buf[i] = bam_nt16_rev_table[buf[i]];
378                 puts((char*)buf);
379                 puts("+");
380                 seq = bam1_qual(b);
381                 for (i = 0; i < qlen; ++i)
382                         buf[i] = 33 + seq[i];
383                 if (b->core.flag & 16) { // reverse
384                         for (i = 0; i < qlen>>1; ++i) {
385                                 int8_t t = buf[qlen - 1 - i];
386                                 buf[qlen - 1 - i] = buf[i];
387                                 buf[i] = t;
388                         }
389                 }
390                 puts((char*)buf);
391         }
392         free(buf);
393         bam_destroy1(b);
394         bam_header_destroy(h);
395         bam_close(fp);
396         return 0;
397 }