]> git.donarmstrong.com Git - samtools.git/blob - bam_sort.c
Merge branch 'master' into mt
[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)
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;
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) fpout = strcmp(out, "-")? bam_open(out, "wu") : bam_dopen(fileno(stdout), "wu");
216         else if (flag & MERGE_LEVEL1) fpout = strcmp(out, "-")? bam_open(out, "w1") : bam_dopen(fileno(stdout), "w1");
217         else fpout = strcmp(out, "-")? bam_open(out, "w") : bam_dopen(fileno(stdout), "w");
218         if (fpout == 0) {
219                 fprintf(stderr, "[%s] fail to create the output file.\n", __func__);
220                 return -1;
221         }
222         bam_header_write(fpout, hout);
223         bam_header_destroy(hout);
224         if (!(flag & MERGE_UNCOMP)) bgzf_mt(fpout, n_threads, 256);
225
226         ks_heapmake(heap, n, heap);
227         while (heap->pos != HEAP_EMPTY) {
228                 bam1_t *b = heap->b;
229                 if (flag & MERGE_RG) {
230                         uint8_t *rg = bam_aux_get(b, "RG");
231                         if (rg) bam_aux_del(b, rg);
232                         bam_aux_append(b, "RG", 'Z', RG_len[heap->i] + 1, (uint8_t*)RG[heap->i]);
233                 }
234                 bam_write1_core(fpout, &b->core, b->data_len, b->data);
235                 if ((j = bam_iter_read(fp[heap->i], iter[heap->i], b)) >= 0) {
236                         heap->pos = ((uint64_t)b->core.tid<<32) | (uint32_t)((int)b->core.pos+1)<<1 | bam1_strand(b);
237                         heap->idx = idx++;
238                 } else if (j == -1) {
239                         heap->pos = HEAP_EMPTY;
240                         free(heap->b->data); free(heap->b);
241                         heap->b = 0;
242                 } else fprintf(stderr, "[bam_merge_core] '%s' is truncated. Continue anyway.\n", fn[heap->i]);
243                 ks_heapadjust(heap, 0, n, heap);
244         }
245
246         if (flag & MERGE_RG) {
247                 for (i = 0; i != n; ++i) free(RG[i]);
248                 free(RG); free(RG_len);
249         }
250         for (i = 0; i != n; ++i) {
251                 bam_iter_destroy(iter[i]);
252                 bam_close(fp[i]);
253         }
254         bam_close(fpout);
255         free(fp); free(heap); free(iter);
256         return 0;
257 }
258
259 int bam_merge_core(int by_qname, const char *out, const char *headers, int n, char * const *fn, int flag, const char *reg)
260 {
261         return bam_merge_core2(by_qname, out, headers, n, fn, flag, reg, 0);
262 }
263
264 int bam_merge(int argc, char *argv[])
265 {
266         int c, is_by_qname = 0, flag = 0, ret = 0, n_threads = 0;
267         char *fn_headers = NULL, *reg = 0;
268
269         while ((c = getopt(argc, argv, "h:nru1R:f@:")) >= 0) {
270                 switch (c) {
271                 case 'r': flag |= MERGE_RG; break;
272                 case 'f': flag |= MERGE_FORCE; break;
273                 case 'h': fn_headers = strdup(optarg); break;
274                 case 'n': is_by_qname = 1; break;
275                 case '1': flag |= MERGE_LEVEL1; break;
276                 case 'u': flag |= MERGE_UNCOMP; break;
277                 case 'R': reg = strdup(optarg); break;
278                 case '@': n_threads = atoi(optarg); break;
279                 }
280         }
281         if (optind + 2 >= argc) {
282                 fprintf(stderr, "\n");
283                 fprintf(stderr, "Usage:   samtools merge [-nr] [-h inh.sam] <out.bam> <in1.bam> <in2.bam> [...]\n\n");
284                 fprintf(stderr, "Options: -n       sort by read names\n");
285                 fprintf(stderr, "         -r       attach RG tag (inferred from file names)\n");
286                 fprintf(stderr, "         -u       uncompressed BAM output\n");
287                 fprintf(stderr, "         -f       overwrite the output BAM if exist\n");
288                 fprintf(stderr, "         -1       compress level 1\n");
289                 fprintf(stderr, "         -@ INT   number of BAM compression threads [0]\n");
290                 fprintf(stderr, "         -R STR   merge file in the specified region STR [all]\n");
291                 fprintf(stderr, "         -h FILE  copy the header in FILE to <out.bam> [in1.bam]\n\n");
292                 fprintf(stderr, "Note: Samtools' merge does not reconstruct the @RG dictionary in the header. Users\n");
293                 fprintf(stderr, "      must provide the correct header with -h, or uses Picard which properly maintains\n");
294                 fprintf(stderr, "      the header dictionary in merging.\n\n");
295                 return 1;
296         }
297         if (!(flag & MERGE_FORCE) && strcmp(argv[optind], "-")) {
298                 FILE *fp = fopen(argv[optind], "rb");
299                 if (fp != NULL) {
300                         fclose(fp);
301                         fprintf(stderr, "[%s] File '%s' exists. Please apply '-f' to overwrite. Abort.\n", __func__, argv[optind]);
302                         return 1;
303                 }
304         }
305         if (bam_merge_core2(is_by_qname, argv[optind], fn_headers, argc - optind - 1, argv + optind + 1, flag, reg, n_threads) < 0) ret = 1;
306         free(reg);
307         free(fn_headers);
308         return ret;
309 }
310
311 /***************
312  * BAM sorting *
313  ***************/
314
315 #include <pthread.h>
316
317 typedef bam1_t *bam1_p;
318
319 static inline int bam1_lt(const bam1_p a, const bam1_p b)
320 {
321         if (g_is_by_qname) {
322                 int t = strnum_cmp(bam1_qname(a), bam1_qname(b));
323                 return (t < 0 || (t == 0 && (a->core.flag&0xc0) < (b->core.flag&0xc0)));
324         } 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)));
325 }
326 KSORT_INIT(sort, bam1_p, bam1_lt)
327
328 typedef struct {
329         size_t buf_len;
330         const char *prefix;
331         bam1_p *buf;
332         const bam_header_t *h;
333         int index;
334 } worker_t;
335
336 static void write_buffer(const char *fn, const char *mode, size_t l, bam1_p *buf, const bam_header_t *h, int n_threads)
337 {
338         size_t i;
339         bamFile fp;
340         fp = strcmp(fn, "-")? bam_open(fn, mode) : bam_dopen(fileno(stdout), mode);
341         if (fp == 0) return;
342         bam_header_write(fp, h);
343         if (n_threads > 1) bgzf_mt(fp, n_threads, 256);
344         for (i = 0; i < l; ++i)
345                 bam_write1_core(fp, &buf[i]->core, buf[i]->data_len, buf[i]->data);
346         bam_close(fp);
347 }
348
349 static void *worker(void *data)
350 {
351         worker_t *w = (worker_t*)data;
352         char *name;
353         ks_mergesort(sort, w->buf_len, w->buf, 0);
354         name = (char*)calloc(strlen(w->prefix) + 20, 1);
355         sprintf(name, "%s.%.4d.bam", w->prefix, w->index);
356         write_buffer(name, "w1", w->buf_len, w->buf, w->h, 0);
357         free(name);
358         return 0;
359 }
360
361 static int sort_blocks(int n_files, size_t k, bam1_p *buf, const char *prefix, const bam_header_t *h, int n_threads)
362 {
363         int i;
364         size_t rest;
365         bam1_p *b;
366         pthread_t *tid;
367         pthread_attr_t attr;
368         worker_t *w;
369
370         if (n_threads < 1) n_threads = 1;
371         if (k < n_threads * 64) n_threads = 1; // use a single thread if we only sort a small batch of records
372         pthread_attr_init(&attr);
373         pthread_attr_setdetachstate(&attr, PTHREAD_CREATE_JOINABLE);
374         w = calloc(n_threads, sizeof(worker_t));
375         tid = calloc(n_threads, sizeof(pthread_t));
376         b = buf; rest = k;
377         for (i = 0; i < n_threads; ++i) {
378                 w[i].buf_len = rest / (n_threads - i);
379                 w[i].buf = b;
380                 w[i].prefix = prefix;
381                 w[i].h = h;
382                 w[i].index = n_files + i;
383                 b += w[i].buf_len; rest -= w[i].buf_len;
384                 pthread_create(&tid[i], &attr, worker, &w[i]);
385         }
386         for (i = 0; i < n_threads; ++i) pthread_join(tid[i], 0);
387         free(tid); free(w);
388         return n_files + n_threads;
389 }
390
391 /*!
392   @abstract Sort an unsorted BAM file based on the chromosome order
393   and the leftmost position of an alignment
394
395   @param  is_by_qname whether to sort by query name
396   @param  fn       name of the file to be sorted
397   @param  prefix   prefix of the output and the temporary files; upon
398                            sucessess, prefix.bam will be written.
399   @param  max_mem  approxiate maximum memory (very inaccurate)
400
401   @discussion It may create multiple temporary subalignment files
402   and then merge them by calling bam_merge_core(). This function is
403   NOT thread safe.
404  */
405 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)
406 {
407         int ret, i, extra_mem, n_files = 0;
408         size_t mem, max_k, k, max_mem;
409         bam_header_t *header;
410         bamFile fp;
411         bam1_t *b, **buf;
412         char *fnout = 0;
413
414         if (n_threads < 2) n_threads = 1;
415         extra_mem = sizeof(void*) + sizeof(void*) + (sizeof(bam1_t) - sizeof(bam1_core_t));
416         g_is_by_qname = is_by_qname;
417         max_k = k = 0; mem = 0;
418         max_mem = _max_mem * n_threads;
419         buf = 0;
420         fp = strcmp(fn, "-")? bam_open(fn, "r") : bam_dopen(fileno(stdin), "r");
421         if (fp == 0) {
422                 fprintf(stderr, "[bam_sort_core] fail to open file %s\n", fn);
423                 return;
424         }
425         header = bam_header_read(fp);
426         // write sub files
427         for (;;) {
428                 if (k == max_k) {
429                         size_t old_max = max_k;
430                         max_k = max_k? max_k<<1 : 0x10000;
431                         buf = realloc(buf, max_k * sizeof(void*));
432                         memset(buf + old_max, 0, sizeof(void*) * (max_k - old_max));
433                 }
434                 if (buf[k] == 0) buf[k] = (bam1_t*)calloc(1, sizeof(bam1_t));
435                 b = buf[k];
436                 if ((ret = bam_read1(fp, b)) < 0) break;
437                 mem += ret + extra_mem;
438                 ++k;
439                 if (mem >= max_mem) {
440                         n_files = sort_blocks(n_files, k, buf, prefix, header, n_threads);
441                         mem = k = 0;
442                 }
443         }
444         if (ret != -1)
445                 fprintf(stderr, "[bam_sort_core] truncated file. Continue anyway.\n");
446         // output file name
447         fnout = calloc(strlen(prefix) + 20, 1);
448         if (is_stdout) sprintf(fnout, "-");
449         else sprintf(fnout, "%s.bam", prefix);
450         // write the final output
451         if (n_files == 0) { // a single block
452                 ks_mergesort(sort, k, buf, 0);
453                 write_buffer(fnout, "w", k, buf, header, n_threads);
454         } else { // then merge
455                 char **fns;
456                 n_files = sort_blocks(n_files, k, buf, prefix, header, n_threads);
457                 fprintf(stderr, "[bam_sort_core] merging from %d files...\n", n_files);
458                 fns = (char**)calloc(n_files, sizeof(char*));
459                 for (i = 0; i < n_files; ++i) {
460                         fns[i] = (char*)calloc(strlen(prefix) + 20, 1);
461                         sprintf(fns[i], "%s.%.4d.bam", prefix, i);
462                 }
463                 bam_merge_core2(is_by_qname, fnout, 0, n_files, fns, 0, 0, n_threads);
464                 for (i = 0; i < n_files; ++i) {
465                         unlink(fns[i]);
466                         free(fns[i]);
467                 }
468                 free(fns);
469         }
470         free(fnout);
471         // free
472         for (k = 0; k < max_k; ++k) {
473                 if (!buf[k]) continue;
474                 free(buf[k]->data);
475                 free(buf[k]);
476         }
477         free(buf);
478         bam_header_destroy(header);
479         bam_close(fp);
480 }
481
482 void bam_sort_core(int is_by_qname, const char *fn, const char *prefix, size_t max_mem)
483 {
484         bam_sort_core_ext(is_by_qname, fn, prefix, max_mem, 0, 0);
485 }
486
487 int bam_sort(int argc, char *argv[])
488 {
489         size_t max_mem = 768<<20; // 512MB
490         int c, is_by_qname = 0, is_stdout = 0, n_threads = 0;
491         while ((c = getopt(argc, argv, "nom:@:")) >= 0) {
492                 switch (c) {
493                 case 'o': is_stdout = 1; break;
494                 case 'n': is_by_qname = 1; break;
495                 case 'm': {
496                                 char *q;
497                                 max_mem = strtol(optarg, &q, 0);
498                                 if (*q == 'k' || *q == 'K') max_mem <<= 10;
499                                 else if (*q == 'm' || *q == 'M') max_mem <<= 20;
500                                 else if (*q == 'g' || *q == 'G') max_mem <<= 30;
501                                 break;
502                         }
503                 case '@': n_threads = atoi(optarg); break;
504                 }
505         }
506         if (optind + 2 > argc) {
507                 fprintf(stderr, "Usage: samtools sort [-on] [-m maxMem=1G] <in.bam> <out.prefix>\n");
508                 return 1;
509         }
510         bam_sort_core_ext(is_by_qname, argv[optind], argv[optind+1], max_mem, is_stdout, n_threads);
511         return 0;
512 }