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