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