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