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