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