New upstream release (no new copyright notices).
[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 inline int strnum_cmp(const char *a, const char *b)
14 {
15         char *pa, *pb;
16         pa = (char*)a; pb = (char*)b;
17         while (*pa && *pb) {
18                 if (isdigit(*pa) && isdigit(*pb)) {
19                         long ai, bi;
20                         ai = strtol(pa, &pa, 10);
21                         bi = strtol(pb, &pb, 10);
22                         if (ai != bi) return ai<bi? -1 : ai>bi? 1 : 0;
23                 } else {
24                         if (*pa != *pb) break;
25                         ++pa; ++pb;
26                 }
27         }
28         if (*pa == *pb)
29                 return (pa-a) < (pb-b)? -1 : (pa-a) > (pb-b)? 1 : 0;
30         return *pa<*pb? -1 : *pa>*pb? 1 : 0;
31 }
32
33 #define HEAP_EMPTY 0xffffffffffffffffull
34
35 typedef struct {
36         int i;
37         uint64_t pos, idx;
38         bam1_t *b;
39 } heap1_t;
40
41 #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))))
42
43 static inline int heap_lt(const heap1_t a, const heap1_t b)
44 {
45         if (g_is_by_qname) {
46                 int t;
47                 if (a.b == 0 || b.b == 0) return a.b == 0? 1 : 0;
48                 t = strnum_cmp(bam1_qname(a.b), bam1_qname(b.b));
49                 return (t > 0 || (t == 0 && __pos_cmp(a, b)));
50         } else return __pos_cmp(a, b);
51 }
52
53 KSORT_INIT(heap, heap1_t, heap_lt)
54
55 static void swap_header_targets(bam_header_t *h1, bam_header_t *h2)
56 {
57         bam_header_t t;
58         t.n_targets = h1->n_targets, h1->n_targets = h2->n_targets, h2->n_targets = t.n_targets;
59         t.target_name = h1->target_name, h1->target_name = h2->target_name, h2->target_name = t.target_name;
60         t.target_len = h1->target_len, h1->target_len = h2->target_len, h2->target_len = t.target_len;
61 }
62
63 static void swap_header_text(bam_header_t *h1, bam_header_t *h2)
64 {
65         int tempi;
66         char *temps;
67         tempi = h1->l_text, h1->l_text = h2->l_text, h2->l_text = tempi;
68         temps = h1->text, h1->text = h2->text, h2->text = temps;
69 }
70
71 #define MERGE_RG     1
72 #define MERGE_UNCOMP 2
73 #define MERGE_LEVEL1 4
74
75 /*!
76   @abstract    Merge multiple sorted BAM.
77   @param  is_by_qname whether to sort by query name
78   @param  out  output BAM file name
79   @param  headers  name of SAM file from which to copy '@' header lines,
80                    or NULL to copy them from the first file to be merged
81   @param  n    number of files to be merged
82   @param  fn   names of files to be merged
83
84   @discussion Padding information may NOT correctly maintained. This
85   function is NOT thread safe.
86  */
87 int bam_merge_core(int by_qname, const char *out, const char *headers, int n, char * const *fn,
88                                         int flag, const char *reg)
89 {
90         bamFile fpout, *fp;
91         heap1_t *heap;
92         bam_header_t *hout = 0;
93         bam_header_t *hheaders = NULL;
94         int i, j, *RG_len = 0;
95         uint64_t idx = 0;
96         char **RG = 0;
97         bam_iter_t *iter = 0;
98
99         if (headers) {
100                 tamFile fpheaders = sam_open(headers);
101                 if (fpheaders == 0) {
102                         const char *message = strerror(errno);
103                         fprintf(stderr, "[bam_merge_core] cannot open '%s': %s\n", headers, message);
104                         return -1;
105                 }
106                 hheaders = sam_header_read(fpheaders);
107                 sam_close(fpheaders);
108         }
109
110         g_is_by_qname = by_qname;
111         fp = (bamFile*)calloc(n, sizeof(bamFile));
112         heap = (heap1_t*)calloc(n, sizeof(heap1_t));
113         iter = (bam_iter_t*)calloc(n, sizeof(bam_iter_t));
114         // prepare RG tag
115         if (flag & MERGE_RG) {
116                 RG = (char**)calloc(n, sizeof(void*));
117                 RG_len = (int*)calloc(n, sizeof(int));
118                 for (i = 0; i != n; ++i) {
119                         int l = strlen(fn[i]);
120                         const char *s = fn[i];
121                         if (l > 4 && strcmp(s + l - 4, ".bam") == 0) l -= 4;
122                         for (j = l - 1; j >= 0; --j) if (s[j] == '/') break;
123                         ++j; l -= j;
124                         RG[i] = calloc(l + 1, 1);
125                         RG_len[i] = l;
126                         strncpy(RG[i], s + j, l);
127                 }
128         }
129         // read the first
130         for (i = 0; i != n; ++i) {
131                 bam_header_t *hin;
132                 fp[i] = bam_open(fn[i], "r");
133                 if (fp[i] == 0) {
134                         int j;
135                         fprintf(stderr, "[bam_merge_core] fail to open file %s\n", fn[i]);
136                         for (j = 0; j < i; ++j) bam_close(fp[j]);
137                         free(fp); free(heap);
138                         // FIXME: possible memory leak
139                         return -1;
140                 }
141                 hin = bam_header_read(fp[i]);
142                 if (i == 0) { // the first BAM
143                         hout = hin;
144                 } else { // validate multiple baf
145                         int min_n_targets = hout->n_targets;
146                         if (hin->n_targets < min_n_targets) min_n_targets = hin->n_targets;
147
148                         for (j = 0; j < min_n_targets; ++j)
149                                 if (strcmp(hout->target_name[j], hin->target_name[j]) != 0) {
150                                         fprintf(stderr, "[bam_merge_core] different target sequence name: '%s' != '%s' in file '%s'\n",
151                                                         hout->target_name[j], hin->target_name[j], fn[i]);
152                                         return -1;
153                                 }
154
155                         // If this input file has additional target reference sequences,
156                         // add them to the headers to be output
157                         if (hin->n_targets > hout->n_targets) {
158                                 swap_header_targets(hout, hin);
159                                 // FIXME Possibly we should also create @SQ text headers
160                                 // for the newly added reference sequences
161                         }
162
163                         bam_header_destroy(hin);
164                 }
165         }
166
167         if (hheaders) {
168                 // If the text headers to be swapped in include any @SQ headers,
169                 // check that they are consistent with the existing binary list
170                 // of reference information.
171                 if (hheaders->n_targets > 0) {
172                         if (hout->n_targets != hheaders->n_targets) {
173                                 fprintf(stderr, "[bam_merge_core] number of @SQ headers in '%s' differs from number of target sequences\n", headers);
174                                 if (!reg) return -1;
175                         }
176                         for (j = 0; j < hout->n_targets; ++j)
177                                 if (strcmp(hout->target_name[j], hheaders->target_name[j]) != 0) {
178                                         fprintf(stderr, "[bam_merge_core] @SQ header '%s' in '%s' differs from target sequence\n", hheaders->target_name[j], headers);
179                                         if (!reg) return -1;
180                                 }
181                 }
182
183                 swap_header_text(hout, hheaders);
184                 bam_header_destroy(hheaders);
185         }
186
187         if (reg) {
188                 int tid, beg, end;
189                 if (bam_parse_region(hout, reg, &tid, &beg, &end) < 0) {
190                         fprintf(stderr, "[%s] Malformated region string or undefined reference name\n", __func__);
191                         return -1;
192                 }
193                 for (i = 0; i < n; ++i) {
194                         bam_index_t *idx;
195                         idx = bam_index_load(fn[i]);
196                         iter[i] = bam_iter_query(idx, tid, beg, end);
197                         bam_index_destroy(idx);
198                 }
199         }
200
201         for (i = 0; i < n; ++i) {
202                 heap1_t *h = heap + i;
203                 h->i = i;
204                 h->b = (bam1_t*)calloc(1, sizeof(bam1_t));
205                 if (bam_iter_read(fp[i], iter[i], h->b) >= 0) {
206                         h->pos = ((uint64_t)h->b->core.tid<<32) | (uint32_t)h->b->core.pos<<1 | bam1_strand(h->b);
207                         h->idx = idx++;
208                 }
209                 else h->pos = HEAP_EMPTY;
210         }
211         if (flag & MERGE_UNCOMP) fpout = strcmp(out, "-")? bam_open(out, "wu") : bam_dopen(fileno(stdout), "wu");
212         else if (flag & MERGE_LEVEL1) fpout = strcmp(out, "-")? bam_open(out, "w1") : bam_dopen(fileno(stdout), "w1");
213         else fpout = strcmp(out, "-")? bam_open(out, "w") : bam_dopen(fileno(stdout), "w");
214         if (fpout == 0) {
215                 fprintf(stderr, "[%s] fail to create the output file.\n", __func__);
216                 return -1;
217         }
218         bam_header_write(fpout, hout);
219         bam_header_destroy(hout);
220
221         ks_heapmake(heap, n, heap);
222         while (heap->pos != HEAP_EMPTY) {
223                 bam1_t *b = heap->b;
224                 if (flag & MERGE_RG) {
225                         uint8_t *rg = bam_aux_get(b, "RG");
226                         if (rg) bam_aux_del(b, rg);
227                         bam_aux_append(b, "RG", 'Z', RG_len[heap->i] + 1, (uint8_t*)RG[heap->i]);
228                 }
229                 bam_write1_core(fpout, &b->core, b->data_len, b->data);
230                 if ((j = bam_iter_read(fp[heap->i], iter[heap->i], b)) >= 0) {
231                         heap->pos = ((uint64_t)b->core.tid<<32) | (uint32_t)b->core.pos<<1 | bam1_strand(b);
232                         heap->idx = idx++;
233                 } else if (j == -1) {
234                         heap->pos = HEAP_EMPTY;
235                         free(heap->b->data); free(heap->b);
236                         heap->b = 0;
237                 } else fprintf(stderr, "[bam_merge_core] '%s' is truncated. Continue anyway.\n", fn[heap->i]);
238                 ks_heapadjust(heap, 0, n, heap);
239         }
240
241         if (flag & MERGE_RG) {
242                 for (i = 0; i != n; ++i) free(RG[i]);
243                 free(RG); free(RG_len);
244         }
245         for (i = 0; i != n; ++i) {
246                 bam_iter_destroy(iter[i]);
247                 bam_close(fp[i]);
248         }
249         bam_close(fpout);
250         free(fp); free(heap); free(iter);
251         return 0;
252 }
253
254 int bam_merge(int argc, char *argv[])
255 {
256         int c, is_by_qname = 0, flag = 0, ret = 0;
257         char *fn_headers = NULL, *reg = 0;
258
259         while ((c = getopt(argc, argv, "h:nru1R:")) >= 0) {
260                 switch (c) {
261                 case 'r': flag |= MERGE_RG; break;
262                 case 'h': fn_headers = strdup(optarg); break;
263                 case 'n': is_by_qname = 1; break;
264                 case '1': flag |= MERGE_LEVEL1; break;
265                 case 'u': flag |= MERGE_UNCOMP; break;
266                 case 'R': reg = strdup(optarg); break;
267                 }
268         }
269         if (optind + 2 >= argc) {
270                 fprintf(stderr, "\n");
271                 fprintf(stderr, "Usage:   samtools merge [-nr] [-h inh.sam] <out.bam> <in1.bam> <in2.bam> [...]\n\n");
272                 fprintf(stderr, "Options: -n       sort by read names\n");
273                 fprintf(stderr, "         -r       attach RG tag (inferred from file names)\n");
274                 fprintf(stderr, "         -u       uncompressed BAM output\n");
275                 fprintf(stderr, "         -1       compress level 1\n");
276                 fprintf(stderr, "         -R STR   merge file in the specified region STR [all]\n");
277                 fprintf(stderr, "         -h FILE  copy the header in FILE to <out.bam> [in1.bam]\n\n");
278                 fprintf(stderr, "Note: Samtools' merge does not reconstruct the @RG dictionary in the header. Users\n");
279                 fprintf(stderr, "      must provide the correct header with -h, or uses Picard which properly maintains\n");
280                 fprintf(stderr, "      the header dictionary in merging.\n\n");
281                 return 1;
282         }
283         if (bam_merge_core(is_by_qname, argv[optind], fn_headers, argc - optind - 1, argv + optind + 1, flag, reg) < 0) ret = 1;
284         free(reg);
285         free(fn_headers);
286         return ret;
287 }
288
289 typedef bam1_t *bam1_p;
290
291 static inline int bam1_lt(const bam1_p a, const bam1_p b)
292 {
293         if (g_is_by_qname) {
294                 int t = strnum_cmp(bam1_qname(a), bam1_qname(b));
295                 return (t < 0 || (t == 0 && (((uint64_t)a->core.tid<<32|a->core.pos) < ((uint64_t)b->core.tid<<32|b->core.pos))));
296         } else return (((uint64_t)a->core.tid<<32|a->core.pos) < ((uint64_t)b->core.tid<<32|b->core.pos));
297 }
298 KSORT_INIT(sort, bam1_p, bam1_lt)
299
300 static void sort_blocks(int n, int k, bam1_p *buf, const char *prefix, const bam_header_t *h, int is_stdout)
301 {
302         char *name, mode[3];
303         int i;
304         bamFile fp;
305         ks_mergesort(sort, k, buf, 0);
306         name = (char*)calloc(strlen(prefix) + 20, 1);
307         if (n >= 0) {
308                 sprintf(name, "%s.%.4d.bam", prefix, n);
309                 strcpy(mode, "w1");
310         } else {
311                 sprintf(name, "%s.bam", prefix);
312                 strcpy(mode, "w");
313         }
314         fp = is_stdout? bam_dopen(fileno(stdout), mode) : bam_open(name, mode);
315         if (fp == 0) {
316                 fprintf(stderr, "[sort_blocks] fail to create file %s.\n", name);
317                 free(name);
318                 // FIXME: possible memory leak
319                 return;
320         }
321         free(name);
322         bam_header_write(fp, h);
323         for (i = 0; i < k; ++i)
324                 bam_write1_core(fp, &buf[i]->core, buf[i]->data_len, buf[i]->data);
325         bam_close(fp);
326 }
327
328 /*!
329   @abstract Sort an unsorted BAM file based on the chromosome order
330   and the leftmost position of an alignment
331
332   @param  is_by_qname whether to sort by query name
333   @param  fn       name of the file to be sorted
334   @param  prefix   prefix of the output and the temporary files; upon
335                            sucessess, prefix.bam will be written.
336   @param  max_mem  approxiate maximum memory (very inaccurate)
337
338   @discussion It may create multiple temporary subalignment files
339   and then merge them by calling bam_merge_core(). This function is
340   NOT thread safe.
341  */
342 void bam_sort_core_ext(int is_by_qname, const char *fn, const char *prefix, size_t max_mem, int is_stdout)
343 {
344         int n, ret, k, i;
345         size_t mem;
346         bam_header_t *header;
347         bamFile fp;
348         bam1_t *b, **buf;
349
350         g_is_by_qname = is_by_qname;
351         n = k = 0; mem = 0;
352         fp = strcmp(fn, "-")? bam_open(fn, "r") : bam_dopen(fileno(stdin), "r");
353         if (fp == 0) {
354                 fprintf(stderr, "[bam_sort_core] fail to open file %s\n", fn);
355                 return;
356         }
357         header = bam_header_read(fp);
358         buf = (bam1_t**)calloc(max_mem / BAM_CORE_SIZE, sizeof(bam1_t*));
359         // write sub files
360         for (;;) {
361                 if (buf[k] == 0) buf[k] = (bam1_t*)calloc(1, sizeof(bam1_t));
362                 b = buf[k];
363                 if ((ret = bam_read1(fp, b)) < 0) break;
364                 mem += ret;
365                 ++k;
366                 if (mem >= max_mem) {
367                         sort_blocks(n++, k, buf, prefix, header, 0);
368                         mem = 0; k = 0;
369                 }
370         }
371         if (ret != -1)
372                 fprintf(stderr, "[bam_sort_core] truncated file. Continue anyway.\n");
373         if (n == 0) sort_blocks(-1, k, buf, prefix, header, is_stdout);
374         else { // then merge
375                 char **fns, *fnout;
376                 fprintf(stderr, "[bam_sort_core] merging from %d files...\n", n+1);
377                 sort_blocks(n++, k, buf, prefix, header, 0);
378                 fnout = (char*)calloc(strlen(prefix) + 20, 1);
379                 if (is_stdout) sprintf(fnout, "-");
380                 else sprintf(fnout, "%s.bam", prefix);
381                 fns = (char**)calloc(n, sizeof(char*));
382                 for (i = 0; i < n; ++i) {
383                         fns[i] = (char*)calloc(strlen(prefix) + 20, 1);
384                         sprintf(fns[i], "%s.%.4d.bam", prefix, i);
385                 }
386                 bam_merge_core(is_by_qname, fnout, 0, n, fns, 0, 0);
387                 free(fnout);
388                 for (i = 0; i < n; ++i) {
389                         unlink(fns[i]);
390                         free(fns[i]);
391                 }
392                 free(fns);
393         }
394         for (k = 0; k < max_mem / BAM_CORE_SIZE; ++k) {
395                 if (buf[k]) {
396                         free(buf[k]->data);
397                         free(buf[k]);
398                 }
399         }
400         free(buf);
401         bam_header_destroy(header);
402         bam_close(fp);
403 }
404
405 void bam_sort_core(int is_by_qname, const char *fn, const char *prefix, size_t max_mem)
406 {
407         bam_sort_core_ext(is_by_qname, fn, prefix, max_mem, 0);
408 }
409
410 int bam_sort(int argc, char *argv[])
411 {
412         size_t max_mem = 500000000;
413         int c, is_by_qname = 0, is_stdout = 0;
414         while ((c = getopt(argc, argv, "nom:")) >= 0) {
415                 switch (c) {
416                 case 'o': is_stdout = 1; break;
417                 case 'n': is_by_qname = 1; break;
418                 case 'm': max_mem = atol(optarg); break;
419                 }
420         }
421         if (optind + 2 > argc) {
422                 fprintf(stderr, "Usage: samtools sort [-on] [-m <maxMem>] <in.bam> <out.prefix>\n");
423                 return 1;
424         }
425         bam_sort_core_ext(is_by_qname, argv[optind], argv[optind+1], max_mem, is_stdout);
426         return 0;
427 }