]> git.donarmstrong.com Git - samtools.git/blob - bam_sort.c
control compression level at command line
[samtools.git] / bam_sort.c
1 #include <stdlib.h>
2 #include <ctype.h>
3 #include <assert.h>
4 #include <errno.h>
5 #include <stdio.h>
6 #include <string.h>
7 #include <unistd.h>
8 #include "bam.h"
9 #include "ksort.h"
10
11 static int g_is_by_qname = 0;
12
13 static int strnum_cmp(const char *_a, const char *_b)
14 {
15         const unsigned char *a = (const unsigned char*)_a, *b = (const unsigned char*)_b;
16         const unsigned char *pa = a, *pb = b;
17         while (*pa && *pb) {
18                 if (isdigit(*pa) && isdigit(*pb)) {
19                         while (*pa == '0') ++pa;
20                         while (*pb == '0') ++pb;
21                         while (isdigit(*pa) && isdigit(*pb) && *pa == *pb) ++pa, ++pb;
22                         if (isdigit(*pa) && isdigit(*pb)) {
23                                 int i = 0;
24                                 while (isdigit(pa[i]) && isdigit(pb[i])) ++i;
25                                 return isdigit(pa[i])? 1 : isdigit(pb[i])? -1 : (int)*pa - (int)*pb;
26                         } else if (isdigit(*pa)) return 1;
27                         else if (isdigit(*pb)) return -1;
28                         else if (pa - a != pb - b) return pa - a < pb - b? 1 : -1;
29                 } else {
30                         if (*pa != *pb) return (int)*pa - (int)*pb;
31                         ++pa; ++pb;
32                 }
33         }
34         return *pa? 1 : *pb? -1 : 0;
35 }
36
37 #define HEAP_EMPTY 0xffffffffffffffffull
38
39 typedef struct {
40         int i;
41         uint64_t pos, idx;
42         bam1_t *b;
43 } heap1_t;
44
45 #define __pos_cmp(a, b) ((a).pos > (b).pos || ((a).pos == (b).pos && ((a).i > (b).i || ((a).i == (b).i && (a).idx > (b).idx))))
46
47 static inline int heap_lt(const heap1_t a, const heap1_t b)
48 {
49         if (g_is_by_qname) {
50                 int t;
51                 if (a.b == 0 || b.b == 0) return a.b == 0? 1 : 0;
52                 t = strnum_cmp(bam1_qname(a.b), bam1_qname(b.b));
53                 return (t > 0 || (t == 0 && (a.b->core.flag&0xc0) > (b.b->core.flag&0xc0)));
54         } else return __pos_cmp(a, b);
55 }
56
57 KSORT_INIT(heap, heap1_t, heap_lt)
58
59 static void swap_header_targets(bam_header_t *h1, bam_header_t *h2)
60 {
61         bam_header_t t;
62         t.n_targets = h1->n_targets, h1->n_targets = h2->n_targets, h2->n_targets = t.n_targets;
63         t.target_name = h1->target_name, h1->target_name = h2->target_name, h2->target_name = t.target_name;
64         t.target_len = h1->target_len, h1->target_len = h2->target_len, h2->target_len = t.target_len;
65 }
66
67 static void swap_header_text(bam_header_t *h1, bam_header_t *h2)
68 {
69         int tempi;
70         char *temps;
71         tempi = h1->l_text, h1->l_text = h2->l_text, h2->l_text = tempi;
72         temps = h1->text, h1->text = h2->text, h2->text = temps;
73 }
74
75 #define MERGE_RG     1
76 #define MERGE_UNCOMP 2
77 #define MERGE_LEVEL1 4
78 #define MERGE_FORCE  8
79
80 /*!
81   @abstract    Merge multiple sorted BAM.
82   @param  is_by_qname whether to sort by query name
83   @param  out  output BAM file name
84   @param  headers  name of SAM file from which to copy '@' header lines,
85                    or NULL to copy them from the first file to be merged
86   @param  n    number of files to be merged
87   @param  fn   names of files to be merged
88
89   @discussion Padding information may NOT correctly maintained. This
90   function is NOT thread safe.
91  */
92 int bam_merge_core2(int by_qname, const char *out, const char *headers, int n, char * const *fn, int flag, const char *reg, int n_threads, int level)
93 {
94         bamFile fpout, *fp;
95         heap1_t *heap;
96         bam_header_t *hout = 0;
97         bam_header_t *hheaders = NULL;
98         int i, j, *RG_len = 0;
99         uint64_t idx = 0;
100         char **RG = 0, mode[8];
101         bam_iter_t *iter = 0;
102
103         if (headers) {
104                 tamFile fpheaders = sam_open(headers);
105                 if (fpheaders == 0) {
106                         const char *message = strerror(errno);
107                         fprintf(stderr, "[bam_merge_core] cannot open '%s': %s\n", headers, message);
108                         return -1;
109                 }
110                 hheaders = sam_header_read(fpheaders);
111                 sam_close(fpheaders);
112         }
113
114         g_is_by_qname = by_qname;
115         fp = (bamFile*)calloc(n, sizeof(bamFile));
116         heap = (heap1_t*)calloc(n, sizeof(heap1_t));
117         iter = (bam_iter_t*)calloc(n, sizeof(bam_iter_t));
118         // prepare RG tag
119         if (flag & MERGE_RG) {
120                 RG = (char**)calloc(n, sizeof(void*));
121                 RG_len = (int*)calloc(n, sizeof(int));
122                 for (i = 0; i != n; ++i) {
123                         int l = strlen(fn[i]);
124                         const char *s = fn[i];
125                         if (l > 4 && strcmp(s + l - 4, ".bam") == 0) l -= 4;
126                         for (j = l - 1; j >= 0; --j) if (s[j] == '/') break;
127                         ++j; l -= j;
128                         RG[i] = calloc(l + 1, 1);
129                         RG_len[i] = l;
130                         strncpy(RG[i], s + j, l);
131                 }
132         }
133         // read the first
134         for (i = 0; i != n; ++i) {
135                 bam_header_t *hin;
136                 fp[i] = bam_open(fn[i], "r");
137                 if (fp[i] == 0) {
138                         int j;
139                         fprintf(stderr, "[bam_merge_core] fail to open file %s\n", fn[i]);
140                         for (j = 0; j < i; ++j) bam_close(fp[j]);
141                         free(fp); free(heap);
142                         // FIXME: possible memory leak
143                         return -1;
144                 }
145                 hin = bam_header_read(fp[i]);
146                 if (i == 0) { // the first BAM
147                         hout = hin;
148                 } else { // validate multiple baf
149                         int min_n_targets = hout->n_targets;
150                         if (hin->n_targets < min_n_targets) min_n_targets = hin->n_targets;
151
152                         for (j = 0; j < min_n_targets; ++j)
153                                 if (strcmp(hout->target_name[j], hin->target_name[j]) != 0) {
154                                         fprintf(stderr, "[bam_merge_core] different target sequence name: '%s' != '%s' in file '%s'\n",
155                                                         hout->target_name[j], hin->target_name[j], fn[i]);
156                                         return -1;
157                                 }
158
159                         // If this input file has additional target reference sequences,
160                         // add them to the headers to be output
161                         if (hin->n_targets > hout->n_targets) {
162                                 swap_header_targets(hout, hin);
163                                 // FIXME Possibly we should also create @SQ text headers
164                                 // for the newly added reference sequences
165                         }
166
167                         bam_header_destroy(hin);
168                 }
169         }
170
171         if (hheaders) {
172                 // If the text headers to be swapped in include any @SQ headers,
173                 // check that they are consistent with the existing binary list
174                 // of reference information.
175                 if (hheaders->n_targets > 0) {
176                         if (hout->n_targets != hheaders->n_targets) {
177                                 fprintf(stderr, "[bam_merge_core] number of @SQ headers in '%s' differs from number of target sequences\n", headers);
178                                 if (!reg) return -1;
179                         }
180                         for (j = 0; j < hout->n_targets; ++j)
181                                 if (strcmp(hout->target_name[j], hheaders->target_name[j]) != 0) {
182                                         fprintf(stderr, "[bam_merge_core] @SQ header '%s' in '%s' differs from target sequence\n", hheaders->target_name[j], headers);
183                                         if (!reg) return -1;
184                                 }
185                 }
186
187                 swap_header_text(hout, hheaders);
188                 bam_header_destroy(hheaders);
189         }
190
191         if (reg) {
192                 int tid, beg, end;
193                 if (bam_parse_region(hout, reg, &tid, &beg, &end) < 0) {
194                         fprintf(stderr, "[%s] Malformated region string or undefined reference name\n", __func__);
195                         return -1;
196                 }
197                 for (i = 0; i < n; ++i) {
198                         bam_index_t *idx;
199                         idx = bam_index_load(fn[i]);
200                         iter[i] = bam_iter_query(idx, tid, beg, end);
201                         bam_index_destroy(idx);
202                 }
203         }
204
205         for (i = 0; i < n; ++i) {
206                 heap1_t *h = heap + i;
207                 h->i = i;
208                 h->b = (bam1_t*)calloc(1, sizeof(bam1_t));
209                 if (bam_iter_read(fp[i], iter[i], h->b) >= 0) {
210                         h->pos = ((uint64_t)h->b->core.tid<<32) | (uint32_t)((int32_t)h->b->core.pos+1)<<1 | bam1_strand(h->b);
211                         h->idx = idx++;
212                 }
213                 else h->pos = HEAP_EMPTY;
214         }
215         if (flag & MERGE_UNCOMP) level = 0;
216         else if (flag & MERGE_LEVEL1) level = 1;
217         strcpy(mode, "w");
218         if (level >= 0) sprintf(mode + 1, "%d", level < 9? level : 9);
219         if ((fpout = strcmp(out, "-")? bam_open(out, "w") : bam_dopen(fileno(stdout), "w")) == 0) {
220                 fprintf(stderr, "[%s] fail to create the output file.\n", __func__);
221                 return -1;
222         }
223         bam_header_write(fpout, hout);
224         bam_header_destroy(hout);
225         if (!(flag & MERGE_UNCOMP)) bgzf_mt(fpout, n_threads, 256);
226
227         ks_heapmake(heap, n, heap);
228         while (heap->pos != HEAP_EMPTY) {
229                 bam1_t *b = heap->b;
230                 if (flag & MERGE_RG) {
231                         uint8_t *rg = bam_aux_get(b, "RG");
232                         if (rg) bam_aux_del(b, rg);
233                         bam_aux_append(b, "RG", 'Z', RG_len[heap->i] + 1, (uint8_t*)RG[heap->i]);
234                 }
235                 bam_write1_core(fpout, &b->core, b->data_len, b->data);
236                 if ((j = bam_iter_read(fp[heap->i], iter[heap->i], b)) >= 0) {
237                         heap->pos = ((uint64_t)b->core.tid<<32) | (uint32_t)((int)b->core.pos+1)<<1 | bam1_strand(b);
238                         heap->idx = idx++;
239                 } else if (j == -1) {
240                         heap->pos = HEAP_EMPTY;
241                         free(heap->b->data); free(heap->b);
242                         heap->b = 0;
243                 } else fprintf(stderr, "[bam_merge_core] '%s' is truncated. Continue anyway.\n", fn[heap->i]);
244                 ks_heapadjust(heap, 0, n, heap);
245         }
246
247         if (flag & MERGE_RG) {
248                 for (i = 0; i != n; ++i) free(RG[i]);
249                 free(RG); free(RG_len);
250         }
251         for (i = 0; i != n; ++i) {
252                 bam_iter_destroy(iter[i]);
253                 bam_close(fp[i]);
254         }
255         bam_close(fpout);
256         free(fp); free(heap); free(iter);
257         return 0;
258 }
259
260 int bam_merge_core(int by_qname, const char *out, const char *headers, int n, char * const *fn, int flag, const char *reg)
261 {
262         return bam_merge_core2(by_qname, out, headers, n, fn, flag, reg, 0, -1);
263 }
264
265 int bam_merge(int argc, char *argv[])
266 {
267         int c, is_by_qname = 0, flag = 0, ret = 0, n_threads = 0, level = -1;
268         char *fn_headers = NULL, *reg = 0;
269
270         while ((c = getopt(argc, argv, "h:nru1R:f@:l:")) >= 0) {
271                 switch (c) {
272                 case 'r': flag |= MERGE_RG; break;
273                 case 'f': flag |= MERGE_FORCE; break;
274                 case 'h': fn_headers = strdup(optarg); break;
275                 case 'n': is_by_qname = 1; break;
276                 case '1': flag |= MERGE_LEVEL1; break;
277                 case 'u': flag |= MERGE_UNCOMP; break;
278                 case 'R': reg = strdup(optarg); break;
279                 case 'l': level = atoi(optarg); break;
280                 case '@': n_threads = atoi(optarg); break;
281                 }
282         }
283         if (optind + 2 >= argc) {
284                 fprintf(stderr, "\n");
285                 fprintf(stderr, "Usage:   samtools merge [-nr] [-h inh.sam] <out.bam> <in1.bam> <in2.bam> [...]\n\n");
286                 fprintf(stderr, "Options: -n       sort by read names\n");
287                 fprintf(stderr, "         -r       attach RG tag (inferred from file names)\n");
288                 fprintf(stderr, "         -u       uncompressed BAM output\n");
289                 fprintf(stderr, "         -f       overwrite the output BAM if exist\n");
290                 fprintf(stderr, "         -1       compress level 1\n");
291                 fprintf(stderr, "         -@ INT   number of BAM compression threads [0]\n");
292                 fprintf(stderr, "         -R STR   merge file in the specified region STR [all]\n");
293                 fprintf(stderr, "         -h FILE  copy the header in FILE to <out.bam> [in1.bam]\n\n");
294                 fprintf(stderr, "Note: Samtools' merge does not reconstruct the @RG dictionary in the header. Users\n");
295                 fprintf(stderr, "      must provide the correct header with -h, or uses Picard which properly maintains\n");
296                 fprintf(stderr, "      the header dictionary in merging.\n\n");
297                 return 1;
298         }
299         if (!(flag & MERGE_FORCE) && strcmp(argv[optind], "-")) {
300                 FILE *fp = fopen(argv[optind], "rb");
301                 if (fp != NULL) {
302                         fclose(fp);
303                         fprintf(stderr, "[%s] File '%s' exists. Please apply '-f' to overwrite. Abort.\n", __func__, argv[optind]);
304                         return 1;
305                 }
306         }
307         if (bam_merge_core2(is_by_qname, argv[optind], fn_headers, argc - optind - 1, argv + optind + 1, flag, reg, n_threads, level) < 0) ret = 1;
308         free(reg);
309         free(fn_headers);
310         return ret;
311 }
312
313 /***************
314  * BAM sorting *
315  ***************/
316
317 #include <pthread.h>
318
319 typedef bam1_t *bam1_p;
320
321 static inline int bam1_lt(const bam1_p a, const bam1_p b)
322 {
323         if (g_is_by_qname) {
324                 int t = strnum_cmp(bam1_qname(a), bam1_qname(b));
325                 return (t < 0 || (t == 0 && (a->core.flag&0xc0) < (b->core.flag&0xc0)));
326         } else return (((uint64_t)a->core.tid<<32|(a->core.pos+1)<<1|bam1_strand(a)) < ((uint64_t)b->core.tid<<32|(b->core.pos+1)<<1|bam1_strand(b)));
327 }
328 KSORT_INIT(sort, bam1_p, bam1_lt)
329
330 typedef struct {
331         size_t buf_len;
332         const char *prefix;
333         bam1_p *buf;
334         const bam_header_t *h;
335         int index;
336 } worker_t;
337
338 static void write_buffer(const char *fn, const char *mode, size_t l, bam1_p *buf, const bam_header_t *h, int n_threads)
339 {
340         size_t i;
341         bamFile fp;
342         fp = strcmp(fn, "-")? bam_open(fn, mode) : bam_dopen(fileno(stdout), mode);
343         if (fp == 0) return;
344         bam_header_write(fp, h);
345         if (n_threads > 1) bgzf_mt(fp, n_threads, 256);
346         for (i = 0; i < l; ++i)
347                 bam_write1_core(fp, &buf[i]->core, buf[i]->data_len, buf[i]->data);
348         bam_close(fp);
349 }
350
351 static void *worker(void *data)
352 {
353         worker_t *w = (worker_t*)data;
354         char *name;
355         ks_mergesort(sort, w->buf_len, w->buf, 0);
356         name = (char*)calloc(strlen(w->prefix) + 20, 1);
357         sprintf(name, "%s.%.4d.bam", w->prefix, w->index);
358         write_buffer(name, "w1", w->buf_len, w->buf, w->h, 0);
359         free(name);
360         return 0;
361 }
362
363 static int sort_blocks(int n_files, size_t k, bam1_p *buf, const char *prefix, const bam_header_t *h, int n_threads)
364 {
365         int i;
366         size_t rest;
367         bam1_p *b;
368         pthread_t *tid;
369         pthread_attr_t attr;
370         worker_t *w;
371
372         if (n_threads < 1) n_threads = 1;
373         if (k < n_threads * 64) n_threads = 1; // use a single thread if we only sort a small batch of records
374         pthread_attr_init(&attr);
375         pthread_attr_setdetachstate(&attr, PTHREAD_CREATE_JOINABLE);
376         w = calloc(n_threads, sizeof(worker_t));
377         tid = calloc(n_threads, sizeof(pthread_t));
378         b = buf; rest = k;
379         for (i = 0; i < n_threads; ++i) {
380                 w[i].buf_len = rest / (n_threads - i);
381                 w[i].buf = b;
382                 w[i].prefix = prefix;
383                 w[i].h = h;
384                 w[i].index = n_files + i;
385                 b += w[i].buf_len; rest -= w[i].buf_len;
386                 pthread_create(&tid[i], &attr, worker, &w[i]);
387         }
388         for (i = 0; i < n_threads; ++i) pthread_join(tid[i], 0);
389         free(tid); free(w);
390         return n_files + n_threads;
391 }
392
393 /*!
394   @abstract Sort an unsorted BAM file based on the chromosome order
395   and the leftmost position of an alignment
396
397   @param  is_by_qname whether to sort by query name
398   @param  fn       name of the file to be sorted
399   @param  prefix   prefix of the output and the temporary files; upon
400                            sucessess, prefix.bam will be written.
401   @param  max_mem  approxiate maximum memory (very inaccurate)
402
403   @discussion It may create multiple temporary subalignment files
404   and then merge them by calling bam_merge_core(). This function is
405   NOT thread safe.
406  */
407 void bam_sort_core_ext(int is_by_qname, const char *fn, const char *prefix, size_t _max_mem, int is_stdout, int n_threads, int level)
408 {
409         int ret, i, extra_mem, n_files = 0;
410         size_t mem, max_k, k, max_mem;
411         bam_header_t *header;
412         bamFile fp;
413         bam1_t *b, **buf;
414         char *fnout = 0;
415
416         if (n_threads < 2) n_threads = 1;
417         extra_mem = sizeof(void*) + sizeof(void*) + (sizeof(bam1_t) - sizeof(bam1_core_t));
418         g_is_by_qname = is_by_qname;
419         max_k = k = 0; mem = 0;
420         max_mem = _max_mem * n_threads;
421         buf = 0;
422         fp = strcmp(fn, "-")? bam_open(fn, "r") : bam_dopen(fileno(stdin), "r");
423         if (fp == 0) {
424                 fprintf(stderr, "[bam_sort_core] fail to open file %s\n", fn);
425                 return;
426         }
427         header = bam_header_read(fp);
428         // write sub files
429         for (;;) {
430                 if (k == max_k) {
431                         size_t old_max = max_k;
432                         max_k = max_k? max_k<<1 : 0x10000;
433                         buf = realloc(buf, max_k * sizeof(void*));
434                         memset(buf + old_max, 0, sizeof(void*) * (max_k - old_max));
435                 }
436                 if (buf[k] == 0) buf[k] = (bam1_t*)calloc(1, sizeof(bam1_t));
437                 b = buf[k];
438                 if ((ret = bam_read1(fp, b)) < 0) break;
439                 mem += ret + extra_mem;
440                 ++k;
441                 if (mem >= max_mem) {
442                         n_files = sort_blocks(n_files, k, buf, prefix, header, n_threads);
443                         mem = k = 0;
444                 }
445         }
446         if (ret != -1)
447                 fprintf(stderr, "[bam_sort_core] truncated file. Continue anyway.\n");
448         // output file name
449         fnout = calloc(strlen(prefix) + 20, 1);
450         if (is_stdout) sprintf(fnout, "-");
451         else sprintf(fnout, "%s.bam", prefix);
452         // write the final output
453         if (n_files == 0) { // a single block
454                 char mode[8];
455                 strcpy(mode, "w");
456                 if (level >= 0) sprintf(mode + 1, "%d", level < 9? level : 9);
457                 ks_mergesort(sort, k, buf, 0);
458                 write_buffer(fnout, mode, k, buf, header, n_threads);
459         } else { // then merge
460                 char **fns;
461                 n_files = sort_blocks(n_files, k, buf, prefix, header, n_threads);
462                 fprintf(stderr, "[bam_sort_core] merging from %d files...\n", n_files);
463                 fns = (char**)calloc(n_files, sizeof(char*));
464                 for (i = 0; i < n_files; ++i) {
465                         fns[i] = (char*)calloc(strlen(prefix) + 20, 1);
466                         sprintf(fns[i], "%s.%.4d.bam", prefix, i);
467                 }
468                 bam_merge_core2(is_by_qname, fnout, 0, n_files, fns, 0, 0, n_threads, level);
469                 for (i = 0; i < n_files; ++i) {
470                         unlink(fns[i]);
471                         free(fns[i]);
472                 }
473                 free(fns);
474         }
475         free(fnout);
476         // free
477         for (k = 0; k < max_k; ++k) {
478                 if (!buf[k]) continue;
479                 free(buf[k]->data);
480                 free(buf[k]);
481         }
482         free(buf);
483         bam_header_destroy(header);
484         bam_close(fp);
485 }
486
487 void bam_sort_core(int is_by_qname, const char *fn, const char *prefix, size_t max_mem)
488 {
489         bam_sort_core_ext(is_by_qname, fn, prefix, max_mem, 0, 0, -1);
490 }
491
492 int bam_sort(int argc, char *argv[])
493 {
494         size_t max_mem = 768<<20; // 512MB
495         int c, is_by_qname = 0, is_stdout = 0, n_threads = 0, level = -1;
496         while ((c = getopt(argc, argv, "nom:@:l:")) >= 0) {
497                 switch (c) {
498                 case 'o': is_stdout = 1; break;
499                 case 'n': is_by_qname = 1; break;
500                 case 'm': {
501                                 char *q;
502                                 max_mem = strtol(optarg, &q, 0);
503                                 if (*q == 'k' || *q == 'K') max_mem <<= 10;
504                                 else if (*q == 'm' || *q == 'M') max_mem <<= 20;
505                                 else if (*q == 'g' || *q == 'G') max_mem <<= 30;
506                                 break;
507                         }
508                 case '@': n_threads = atoi(optarg); break;
509                 case 'l': level = atoi(optarg); break;
510                 }
511         }
512         if (optind + 2 > argc) {
513                 fprintf(stderr, "Usage: samtools sort [-on] [-m maxMem=1G] <in.bam> <out.prefix>\n");
514                 return 1;
515         }
516         bam_sort_core_ext(is_by_qname, argv[optind], argv[optind+1], max_mem, is_stdout, n_threads, level);
517         return 0;
518 }