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