Imported Upstream version 0.2.0
[tabix.git] / index.c
1 #include <ctype.h>
2 #include <assert.h>
3 #include <sys/stat.h>
4 #include "khash.h"
5 #include "ksort.h"
6 #include "kstring.h"
7 #include "bam_endian.h"
8 #ifdef _USE_KNETFILE
9 #include "knetfile.h"
10 #endif
11 #include "tabix.h"
12
13 #define TAD_MIN_CHUNK_GAP 32768
14 // 1<<14 is the size of minimum bin.
15 #define TAD_LIDX_SHIFT    14
16
17 typedef struct {
18         uint64_t u, v;
19 } pair64_t;
20
21 #define pair64_lt(a,b) ((a).u < (b).u)
22 KSORT_INIT(off, pair64_t, pair64_lt)
23
24 typedef struct {
25         uint32_t m, n;
26         pair64_t *list;
27 } ti_binlist_t;
28
29 typedef struct {
30         int32_t n, m;
31         uint64_t *offset;
32 } ti_lidx_t;
33
34 KHASH_MAP_INIT_INT(i, ti_binlist_t)
35 KHASH_MAP_INIT_STR(s, int)
36
37 struct __ti_index_t {
38         ti_conf_t conf;
39         int32_t n, max;
40         khash_t(s) *tname;
41         khash_t(i) **index;
42         ti_lidx_t *index2;
43 };
44
45 struct __ti_iter_t {
46         int from_first; // read from the first record; no random access
47         int tid, beg, end, n_off, i, finished;
48         uint64_t curr_off;
49         kstring_t str;
50         const ti_index_t *idx;
51         pair64_t *off;
52 };
53
54 typedef struct {
55         int tid, beg, end, bin;
56 } ti_intv_t;
57
58 ti_conf_t ti_conf_gff = { 0, 1, 4, 5, '#', 0 };
59 ti_conf_t ti_conf_bed = { TI_FLAG_UCSC, 1, 2, 3, '#', 0 };
60 ti_conf_t ti_conf_psltbl = { TI_FLAG_UCSC, 15, 17, 18, '#', 0 };
61 ti_conf_t ti_conf_sam = { TI_PRESET_SAM, 3, 4, 0, '@', 0 };
62 ti_conf_t ti_conf_vcf = { TI_PRESET_VCF, 1, 2, 0, '#', 0 };
63
64 /***************
65  * read a line *
66  ***************/
67
68 /*
69 int ti_readline(BGZF *fp, kstring_t *str)
70 {
71         int c, l = 0;
72         str->l = 0;
73         while ((c = bgzf_getc(fp)) >= 0 && c != '\n') {
74                 ++l;
75                 if (c != '\r') kputc(c, str);
76         }
77         if (c < 0 && l == 0) return -1; // end of file
78         return str->l;
79 }
80 */
81
82 /* Below is a faster implementation largely equivalent to the one
83  * commented out above. */
84 int ti_readline(BGZF *fp, kstring_t *str)
85 {
86         int l, state = 0;
87         unsigned char *buf = (unsigned char*)fp->uncompressed_block;
88         str->l = 0;
89         do {
90                 if (fp->block_offset >= fp->block_length) {
91                         if (bgzf_read_block(fp) != 0) { state = -2; break; }
92                         if (fp->block_length == 0) { state = -1; break; }
93                 }
94                 for (l = fp->block_offset; l < fp->block_length && buf[l] != '\n'; ++l);
95                 if (l < fp->block_length) state = 1;
96                 l -= fp->block_offset;
97                 if (str->l + l + 1 >= str->m) {
98                         str->m = str->l + l + 2;
99                         kroundup32(str->m);
100                         str->s = (char*)realloc(str->s, str->m);
101                 }
102                 memcpy(str->s + str->l, buf + fp->block_offset, l);
103                 str->l += l;
104                 fp->block_offset += l + 1;
105                 if (fp->block_offset >= fp->block_length) {
106 #ifdef _USE_KNETFILE
107                         fp->block_address = knet_tell(fp->x.fpr);
108 #else
109                         fp->block_address = ftello(fp->file);
110 #endif
111                         fp->block_offset = 0;
112                         fp->block_length = 0;
113                 } 
114         } while (state == 0);
115         if (str->l == 0 && state < 0) return state;
116         str->s[str->l] = 0;
117         return str->l;
118 }
119
120 /*************************************
121  * get the interval from a data line *
122  *************************************/
123
124 static inline int ti_reg2bin(uint32_t beg, uint32_t end)
125 {
126         --end;
127         if (beg>>14 == end>>14) return 4681 + (beg>>14);
128         if (beg>>17 == end>>17) return  585 + (beg>>17);
129         if (beg>>20 == end>>20) return   73 + (beg>>20);
130         if (beg>>23 == end>>23) return    9 + (beg>>23);
131         if (beg>>26 == end>>26) return    1 + (beg>>26);
132         return 0;
133 }
134
135 static int get_tid(ti_index_t *idx, const char *ss)
136 {
137         khint_t k;
138         int tid;
139         k = kh_get(s, idx->tname, ss);
140         if (k == kh_end(idx->tname)) { // a new target sequence
141                 int ret, size;
142                 // update idx->n, ->max, ->index and ->index2
143                 if (idx->n == idx->max) {
144                         idx->max = idx->max? idx->max<<1 : 8;
145                         idx->index = realloc(idx->index, idx->max * sizeof(void*));
146                         idx->index2 = realloc(idx->index2, idx->max * sizeof(ti_lidx_t));
147                 }
148                 memset(&idx->index2[idx->n], 0, sizeof(ti_lidx_t));
149                 idx->index[idx->n++] = kh_init(i);
150                 // update ->tname
151                 tid = size = kh_size(idx->tname);
152                 k = kh_put(s, idx->tname, strdup(ss), &ret);
153                 kh_value(idx->tname, k) = size;
154                 assert(idx->n == kh_size(idx->tname));
155         } else tid = kh_value(idx->tname, k);
156         return tid;
157 }
158
159 static int get_intv(ti_index_t *idx, kstring_t *str, ti_intv_t *intv)
160 {
161         int i, b = 0, id = 1;
162         char *s;
163         intv->tid = intv->beg = intv->end = intv->bin = -1;
164         for (i = 0; i <= str->l; ++i) {
165                 if (str->s[i] == '\t' || str->s[i] == 0) {
166                         if (id == idx->conf.sc) {
167                                 str->s[i] = 0;
168                                 intv->tid = get_tid(idx, str->s + b);
169                                 if (i != str->l) str->s[i] = '\t';
170                         } else if (id == idx->conf.bc) {
171                                 // here ->beg is 0-based.
172                                 intv->beg = intv->end = strtol(str->s + b, &s, 0);
173                                 if (!(idx->conf.preset&TI_FLAG_UCSC)) --intv->beg;
174                                 else ++intv->end;
175                                 if (intv->beg < 0) intv->beg = 0;
176                                 if (intv->end < 1) intv->end = 1;
177                         } else {
178                                 if ((idx->conf.preset&0xffff) == TI_PRESET_GENERIC) {
179                                         if (id == idx->conf.ec) intv->end = strtol(str->s + b, &s, 0);
180                                 } else if ((idx->conf.preset&0xffff) == TI_PRESET_SAM) {
181                                         if (id == 6) { // CIGAR
182                                                 int l = 0, op;
183                                                 char *t;
184                                                 for (s = str->s + b; s < str->s + i;) {
185                                                         long x = strtol(s, &t, 10);
186                                                         op = toupper(*t);
187                                                         if (op == 'M' || op == 'D' || op == 'N') l += x;
188                                                         s = t + 1;
189                                                 }
190                                                 if (l == 0) l = 1;
191                                                 intv->end = intv->beg + l;
192                                         }
193                                 } else if ((idx->conf.preset&0xffff) == TI_PRESET_VCF) {
194                                         // FIXME: the following is NOT tested and is likely to be buggy
195                                         if (id == 5) { // ALT
196                                                 char *t;
197                                                 int max = 1;
198                                                 for (s = str->s + b; s < str->s + i;) {
199                                                         if (s[i] == 'D') {
200                                                                 long x = strtol(s + 1, &t, 10);
201                                                                 if (x > max) max = x;
202                                                                 s = t + 1;
203                                                         } else ++s;
204                                                 }
205                                                 intv->end = intv->beg + max;
206                                         }
207                                 }
208                         }
209                         b = i + 1;
210                         ++id;
211                 }
212         }
213         if (intv->tid < 0 || intv->beg < 0 || intv->end < 0) return -1;
214         intv->bin = ti_reg2bin(intv->beg, intv->end);
215         return 0;
216 }
217
218 /************
219  * indexing *
220  ************/
221
222 // requirement: len <= LEN_MASK
223 static inline void insert_offset(khash_t(i) *h, int bin, uint64_t beg, uint64_t end)
224 {
225         khint_t k;
226         ti_binlist_t *l;
227         int ret;
228         k = kh_put(i, h, bin, &ret);
229         l = &kh_value(h, k);
230         if (ret) { // not present
231                 l->m = 1; l->n = 0;
232                 l->list = (pair64_t*)calloc(l->m, 16);
233         }
234         if (l->n == l->m) {
235                 l->m <<= 1;
236                 l->list = (pair64_t*)realloc(l->list, l->m * 16);
237         }
238         l->list[l->n].u = beg; l->list[l->n++].v = end;
239 }
240
241 static inline void insert_offset2(ti_lidx_t *index2, int _beg, int _end, uint64_t offset)
242 {
243         int i, beg, end;
244         beg = _beg >> TAD_LIDX_SHIFT;
245         end = (_end - 1) >> TAD_LIDX_SHIFT;
246         if (index2->m < end + 1) {
247                 int old_m = index2->m;
248                 index2->m = end + 1;
249                 kroundup32(index2->m);
250                 index2->offset = (uint64_t*)realloc(index2->offset, index2->m * 8);
251                 memset(index2->offset + old_m, 0, 8 * (index2->m - old_m));
252         }
253         if (beg == end) {
254                 if (index2->offset[beg] == 0) index2->offset[beg] = offset;
255         } else {
256                 for (i = beg; i <= end; ++i)
257                         if (index2->offset[i] == 0) index2->offset[i] = offset;
258         }
259         if (index2->n < end + 1) index2->n = end + 1;
260 }
261
262 static void merge_chunks(ti_index_t *idx)
263 {
264         khash_t(i) *index;
265         int i, l, m;
266         khint_t k;
267         for (i = 0; i < idx->n; ++i) {
268                 index = idx->index[i];
269                 for (k = kh_begin(index); k != kh_end(index); ++k) {
270                         ti_binlist_t *p;
271                         if (!kh_exist(index, k)) continue;
272                         p = &kh_value(index, k);
273                         m = 0;
274                         for (l = 1; l < p->n; ++l) {
275                                 if (p->list[m].v>>16 == p->list[l].u>>16) p->list[m].v = p->list[l].v;
276                                 else p->list[++m] = p->list[l];
277                         } // ~for(l)
278                         p->n = m + 1;
279                 } // ~for(k)
280         } // ~for(i)
281 }
282
283 static void fill_missing(ti_index_t *idx)
284 {
285         int i, j;
286         for (i = 0; i < idx->n; ++i) {
287                 ti_lidx_t *idx2 = &idx->index2[i];
288                 for (j = 1; j < idx2->n; ++j)
289                         if (idx2->offset[j] == 0)
290                                 idx2->offset[j] = idx2->offset[j-1];
291         }
292 }
293
294 ti_index_t *ti_index_core(BGZF *fp, const ti_conf_t *conf)
295 {
296         int ret;
297         ti_index_t *idx;
298         uint32_t last_bin, save_bin;
299         int32_t last_coor, last_tid, save_tid;
300         uint64_t save_off, last_off, lineno = 0;
301         kstring_t *str;
302
303         str = calloc(1, sizeof(kstring_t));
304
305         idx = (ti_index_t*)calloc(1, sizeof(ti_index_t));
306         idx->conf = *conf;
307         idx->n = idx->max = 0;
308         idx->tname = kh_init(s);
309         idx->index = 0;
310         idx->index2 = 0;
311
312         save_bin = save_tid = last_tid = last_bin = 0xffffffffu;
313         save_off = last_off = bgzf_tell(fp); last_coor = 0xffffffffu;
314         while ((ret = ti_readline(fp, str)) >= 0) {
315                 ti_intv_t intv;
316                 ++lineno;
317                 if (lineno <= idx->conf.line_skip || str->s[0] == idx->conf.meta_char) {
318                         last_off = bgzf_tell(fp);
319                         continue;
320                 }
321                 get_intv(idx, str, &intv);
322                 if (last_tid != intv.tid) { // change of chromosomes
323                         last_tid = intv.tid;
324                         last_bin = 0xffffffffu;
325                 } else if (last_coor > intv.beg) {
326                         fprintf(stderr, "[ti_index_core] the file out of order at line %llu\n", (unsigned long long)lineno);
327                         exit(1);
328                 }
329                 insert_offset2(&idx->index2[intv.tid], intv.beg, intv.end, last_off);
330                 if (intv.bin != last_bin) { // then possibly write the binning index
331                         if (save_bin != 0xffffffffu) // save_bin==0xffffffffu only happens to the first record
332                                 insert_offset(idx->index[save_tid], save_bin, save_off, last_off);
333                         save_off = last_off;
334                         save_bin = last_bin = intv.bin;
335                         save_tid = intv.tid;
336                         if (save_tid < 0) break;
337                 }
338                 if (bgzf_tell(fp) <= last_off) {
339                         fprintf(stderr, "[ti_index_core] bug in BGZF: %llx < %llx\n",
340                                         (unsigned long long)bgzf_tell(fp), (unsigned long long)last_off);
341                         exit(1);
342                 }
343                 last_off = bgzf_tell(fp);
344                 last_coor = intv.beg;
345         }
346         if (save_tid >= 0) insert_offset(idx->index[save_tid], save_bin, save_off, bgzf_tell(fp));
347         merge_chunks(idx);
348         fill_missing(idx);
349
350         free(str->s); free(str);
351         return idx;
352 }
353
354 void ti_index_destroy(ti_index_t *idx)
355 {
356         khint_t k;
357         int i;
358         if (idx == 0) return;
359         // destroy the name hash table
360         for (k = kh_begin(idx->tname); k != kh_end(idx->tname); ++k) {
361                 if (kh_exist(idx->tname, k))
362                         free((char*)kh_key(idx->tname, k));
363         }
364         kh_destroy(s, idx->tname);
365         // destroy the binning index
366         for (i = 0; i < idx->n; ++i) {
367                 khash_t(i) *index = idx->index[i];
368                 ti_lidx_t *index2 = idx->index2 + i;
369                 for (k = kh_begin(index); k != kh_end(index); ++k) {
370                         if (kh_exist(index, k))
371                                 free(kh_value(index, k).list);
372                 }
373                 kh_destroy(i, index);
374                 free(index2->offset);
375         }
376         free(idx->index);
377         // destroy the linear index
378         free(idx->index2);
379         free(idx);
380 }
381
382 /******************
383  * index file I/O *
384  ******************/
385
386 void ti_index_save(const ti_index_t *idx, BGZF *fp)
387 {
388         int32_t i, size, ti_is_be;
389         khint_t k;
390         ti_is_be = bam_is_big_endian();
391         bgzf_write(fp, "TBI\1", 4);
392         if (ti_is_be) {
393                 uint32_t x = idx->n;
394                 bgzf_write(fp, bam_swap_endian_4p(&x), 4);
395         } else bgzf_write(fp, &idx->n, 4);
396         assert(sizeof(ti_conf_t) == 24);
397         if (ti_is_be) { // write ti_conf_t;
398                 uint32_t x[6];
399                 memcpy(x, &idx->conf, 24);
400                 for (i = 0; i < 6; ++i) bgzf_write(fp, bam_swap_endian_4p(&x[i]), 4);
401         } else bgzf_write(fp, &idx->conf, sizeof(ti_conf_t));
402         { // write target names
403                 char **name;
404                 int32_t l = 0;
405                 name = calloc(kh_size(idx->tname), sizeof(void*));
406                 for (k = kh_begin(idx->tname); k != kh_end(idx->tname); ++k)
407                         if (kh_exist(idx->tname, k))
408                                 name[kh_value(idx->tname, k)] = (char*)kh_key(idx->tname, k);
409                 for (i = 0; i < kh_size(idx->tname); ++i)
410                         l += strlen(name[i]) + 1;
411                 if (ti_is_be) bgzf_write(fp, bam_swap_endian_4p(&l), 4);
412                 else bgzf_write(fp, &l, 4);
413                 for (i = 0; i < kh_size(idx->tname); ++i)
414                         bgzf_write(fp, name[i], strlen(name[i]) + 1);
415                 free(name);
416         }
417         for (i = 0; i < idx->n; ++i) {
418                 khash_t(i) *index = idx->index[i];
419                 ti_lidx_t *index2 = idx->index2 + i;
420                 // write binning index
421                 size = kh_size(index);
422                 if (ti_is_be) { // big endian
423                         uint32_t x = size;
424                         bgzf_write(fp, bam_swap_endian_4p(&x), 4);
425                 } else bgzf_write(fp, &size, 4);
426                 for (k = kh_begin(index); k != kh_end(index); ++k) {
427                         if (kh_exist(index, k)) {
428                                 ti_binlist_t *p = &kh_value(index, k);
429                                 if (ti_is_be) { // big endian
430                                         uint32_t x;
431                                         x = kh_key(index, k); bgzf_write(fp, bam_swap_endian_4p(&x), 4);
432                                         x = p->n; bgzf_write(fp, bam_swap_endian_4p(&x), 4);
433                                         for (x = 0; (int)x < p->n; ++x) {
434                                                 bam_swap_endian_8p(&p->list[x].u);
435                                                 bam_swap_endian_8p(&p->list[x].v);
436                                         }
437                                         bgzf_write(fp, p->list, 16 * p->n);
438                                         for (x = 0; (int)x < p->n; ++x) {
439                                                 bam_swap_endian_8p(&p->list[x].u);
440                                                 bam_swap_endian_8p(&p->list[x].v);
441                                         }
442                                 } else {
443                                         bgzf_write(fp, &kh_key(index, k), 4);
444                                         bgzf_write(fp, &p->n, 4);
445                                         bgzf_write(fp, p->list, 16 * p->n);
446                                 }
447                         }
448                 }
449                 // write linear index (index2)
450                 if (ti_is_be) {
451                         int x = index2->n;
452                         bgzf_write(fp, bam_swap_endian_4p(&x), 4);
453                 } else bgzf_write(fp, &index2->n, 4);
454                 if (ti_is_be) { // big endian
455                         int x;
456                         for (x = 0; (int)x < index2->n; ++x)
457                                 bam_swap_endian_8p(&index2->offset[x]);
458                         bgzf_write(fp, index2->offset, 8 * index2->n);
459                         for (x = 0; (int)x < index2->n; ++x)
460                                 bam_swap_endian_8p(&index2->offset[x]);
461                 } else bgzf_write(fp, index2->offset, 8 * index2->n);
462         }
463 }
464
465 static ti_index_t *ti_index_load_core(BGZF *fp)
466 {
467         int i, ti_is_be;
468         char magic[4];
469         ti_index_t *idx;
470         ti_is_be = bam_is_big_endian();
471         if (fp == 0) {
472                 fprintf(stderr, "[ti_index_load_core] fail to load index.\n");
473                 return 0;
474         }
475         bgzf_read(fp, magic, 4);
476         if (strncmp(magic, "TBI\1", 4)) {
477                 fprintf(stderr, "[ti_index_load] wrong magic number.\n");
478                 return 0;
479         }
480         idx = (ti_index_t*)calloc(1, sizeof(ti_index_t));       
481         bgzf_read(fp, &idx->n, 4);
482         if (ti_is_be) bam_swap_endian_4p(&idx->n);
483         idx->tname = kh_init(s);
484         idx->index = (khash_t(i)**)calloc(idx->n, sizeof(void*));
485         idx->index2 = (ti_lidx_t*)calloc(idx->n, sizeof(ti_lidx_t));
486         // read idx->conf
487         bgzf_read(fp, &idx->conf, sizeof(ti_conf_t));
488         if (ti_is_be) {
489                 bam_swap_endian_4p(&idx->conf.preset);
490                 bam_swap_endian_4p(&idx->conf.sc);
491                 bam_swap_endian_4p(&idx->conf.bc);
492                 bam_swap_endian_4p(&idx->conf.ec);
493                 bam_swap_endian_4p(&idx->conf.meta_char);
494                 bam_swap_endian_4p(&idx->conf.line_skip);
495         }
496         { // read target names
497                 int j, ret;
498                 kstring_t *str;
499                 int32_t l;
500                 uint8_t *buf;
501                 bgzf_read(fp, &l, 4);
502                 if (ti_is_be) bam_swap_endian_4p(&l);
503                 buf = calloc(l, 1);
504                 bgzf_read(fp, buf, l);
505                 str = calloc(1, sizeof(kstring_t));
506                 for (i = j = 0; i < l; ++i) {
507                         if (buf[i] == 0) {
508                                 khint_t k = kh_put(s, idx->tname, strdup(str->s), &ret);
509                                 kh_value(idx->tname, k) = j++;
510                                 str->l = 0;
511                         } else kputc(buf[i], str);
512                 }
513                 free(str->s); free(str); free(buf);
514         }
515         for (i = 0; i < idx->n; ++i) {
516                 khash_t(i) *index;
517                 ti_lidx_t *index2 = idx->index2 + i;
518                 uint32_t key, size;
519                 khint_t k;
520                 int j, ret;
521                 ti_binlist_t *p;
522                 index = idx->index[i] = kh_init(i);
523                 // load binning index
524                 bgzf_read(fp, &size, 4);
525                 if (ti_is_be) bam_swap_endian_4p(&size);
526                 for (j = 0; j < (int)size; ++j) {
527                         bgzf_read(fp, &key, 4);
528                         if (ti_is_be) bam_swap_endian_4p(&key);
529                         k = kh_put(i, index, key, &ret);
530                         p = &kh_value(index, k);
531                         bgzf_read(fp, &p->n, 4);
532                         if (ti_is_be) bam_swap_endian_4p(&p->n);
533                         p->m = p->n;
534                         p->list = (pair64_t*)malloc(p->m * 16);
535                         bgzf_read(fp, p->list, 16 * p->n);
536                         if (ti_is_be) {
537                                 int x;
538                                 for (x = 0; x < p->n; ++x) {
539                                         bam_swap_endian_8p(&p->list[x].u);
540                                         bam_swap_endian_8p(&p->list[x].v);
541                                 }
542                         }
543                 }
544                 // load linear index
545                 bgzf_read(fp, &index2->n, 4);
546                 if (ti_is_be) bam_swap_endian_4p(&index2->n);
547                 index2->m = index2->n;
548                 index2->offset = (uint64_t*)calloc(index2->m, 8);
549                 bgzf_read(fp, index2->offset, index2->n * 8);
550                 if (ti_is_be)
551                         for (j = 0; j < index2->n; ++j) bam_swap_endian_8p(&index2->offset[j]);
552         }
553         return idx;
554 }
555
556 ti_index_t *ti_index_load_local(const char *fnidx)
557 {
558         BGZF *fp;
559         fp = bgzf_open(fnidx, "r");
560         if (fp) {
561                 ti_index_t *idx = ti_index_load_core(fp);
562                 bgzf_close(fp);
563                 return idx;
564         } else return 0;
565 }
566
567 #ifdef _USE_KNETFILE
568 static void download_from_remote(const char *url)
569 {
570         const int buf_size = 1 * 1024 * 1024;
571         char *fn;
572         FILE *fp;
573         uint8_t *buf;
574         knetFile *fp_remote;
575         int l;
576         if (strstr(url, "ftp://") != url && strstr(url, "http://") != url) return;
577         l = strlen(url);
578         for (fn = (char*)url + l - 1; fn >= url; --fn)
579                 if (*fn == '/') break;
580         ++fn; // fn now points to the file name
581         fp_remote = knet_open(url, "r");
582         if (fp_remote == 0) {
583                 fprintf(stderr, "[download_from_remote] fail to open remote file.\n");
584                 return;
585         }
586         if ((fp = fopen(fn, "w")) == 0) {
587                 fprintf(stderr, "[download_from_remote] fail to create file in the working directory.\n");
588                 knet_close(fp_remote);
589                 return;
590         }
591         buf = (uint8_t*)calloc(buf_size, 1);
592         while ((l = knet_read(fp_remote, buf, buf_size)) != 0)
593                 fwrite(buf, 1, l, fp);
594         free(buf);
595         fclose(fp);
596         knet_close(fp_remote);
597 }
598 #else
599 static void download_from_remote(const char *url)
600 {
601         return;
602 }
603 #endif
604
605 static char *get_local_version(const char *fn)
606 {
607     struct stat sbuf;
608         char *fnidx = (char*)calloc(strlen(fn) + 5, 1);
609         strcat(strcpy(fnidx, fn), ".tbi");
610         if ((strstr(fnidx, "ftp://") == fnidx || strstr(fnidx, "http://") == fnidx)) {
611                 char *p, *url;
612                 int l = strlen(fnidx);
613                 for (p = fnidx + l - 1; p >= fnidx; --p)
614                         if (*p == '/') break;
615                 url = fnidx; fnidx = strdup(p + 1);
616                 if (stat(fnidx, &sbuf) == 0) {
617                         free(url);
618                         return fnidx;
619                 }
620                 fprintf(stderr, "[%s] downloading the index file...\n", __func__);
621                 download_from_remote(url);
622                 free(url);
623         }
624     if (stat(fnidx, &sbuf) == 0) return fnidx;
625         free(fnidx); return 0;
626 }
627
628 const char **ti_seqname(const ti_index_t *idx, int *n)
629 {
630         const char **names;
631         khint_t k;
632         *n = idx->n;
633         names = calloc(idx->n, sizeof(void*));
634         for (k = kh_begin(idx->tname); k < kh_end(idx->tname); ++k)
635                 if (kh_exist(idx->tname, k))
636                         names[kh_val(idx->tname, k)] = kh_key(idx->tname, k);
637         return names;
638 }
639
640 ti_index_t *ti_index_load(const char *fn)
641 {
642         ti_index_t *idx;
643     char *fname = get_local_version(fn);
644         if (fname == 0) return 0;
645         idx = ti_index_load_local(fname);
646     free(fname);
647         if (idx == 0) fprintf(stderr, "[ti_index_load] fail to load BAM index.\n");
648         return idx;
649 }
650
651 int ti_index_build2(const char *fn, const ti_conf_t *conf, const char *_fnidx)
652 {
653         char *fnidx;
654         BGZF *fp, *fpidx;
655         ti_index_t *idx;
656         if ((fp = bgzf_open(fn, "r")) == 0) {
657                 fprintf(stderr, "[ti_index_build2] fail to open the BAM file.\n");
658                 return -1;
659         }
660         idx = ti_index_core(fp, conf);
661         bgzf_close(fp);
662         if (_fnidx == 0) {
663                 fnidx = (char*)calloc(strlen(fn) + 5, 1);
664                 strcpy(fnidx, fn); strcat(fnidx, ".tbi");
665         } else fnidx = strdup(_fnidx);
666         fpidx = bgzf_open(fnidx, "w");
667         if (fpidx == 0) {
668                 fprintf(stderr, "[ti_index_build2] fail to create the index file.\n");
669                 free(fnidx);
670                 return -1;
671         }
672         ti_index_save(idx, fpidx);
673         ti_index_destroy(idx);
674         bgzf_close(fpidx);
675         free(fnidx);
676         return 0;
677 }
678
679 int ti_index_build(const char *fn, const ti_conf_t *conf)
680 {
681         return ti_index_build2(fn, conf, 0);
682 }
683
684 /********************************************
685  * parse a region in the format chr:beg-end *
686  ********************************************/
687
688 int ti_get_tid(const ti_index_t *idx, const char *name)
689 {
690         khiter_t iter;
691         const khash_t(s) *h = idx->tname;
692         iter = kh_get(s, h, name); /* get the tid */
693         if (iter == kh_end(h)) return -1;
694         return kh_value(h, iter);
695 }
696
697 int ti_parse_region(const ti_index_t *idx, const char *str, int *tid, int *begin, int *end)
698 {
699         char *s, *p;
700         int i, l, k;
701         l = strlen(str);
702         p = s = (char*)malloc(l+1);
703         /* squeeze out "," */
704         for (i = k = 0; i != l; ++i)
705                 if (str[i] != ',' && !isspace(str[i])) s[k++] = str[i];
706         s[k] = 0;
707         for (i = 0; i != k; ++i) if (s[i] == ':') break;
708         s[i] = 0;
709         if ((*tid = ti_get_tid(idx, s)) < 0) {
710                 free(s);
711                 return -1;
712         }
713         if (i == k) { /* dump the whole sequence */
714                 *begin = 0; *end = 1<<29; free(s);
715                 return 0;
716         }
717         for (p = s + i + 1; i != k; ++i) if (s[i] == '-') break;
718         *begin = atoi(p);
719         if (i < k) {
720                 p = s + i + 1;
721                 *end = atoi(p);
722         } else *end = 1<<29;
723         if (*begin > 0) --*begin;
724         free(s);
725         if (*begin > *end) return -1;
726         return 0;
727 }
728
729 /*******************************
730  * retrieve a specified region *
731  *******************************/
732
733 #define MAX_BIN 37450 // =(8^6-1)/7+1
734
735 static inline int reg2bins(uint32_t beg, uint32_t end, uint16_t list[MAX_BIN])
736 {
737         int i = 0, k;
738         if (beg >= end) return 0;
739         if (end >= 1u<<29) end = 1u<<29;
740         --end;
741         list[i++] = 0;
742         for (k =    1 + (beg>>26); k <=    1 + (end>>26); ++k) list[i++] = k;
743         for (k =    9 + (beg>>23); k <=    9 + (end>>23); ++k) list[i++] = k;
744         for (k =   73 + (beg>>20); k <=   73 + (end>>20); ++k) list[i++] = k;
745         for (k =  585 + (beg>>17); k <=  585 + (end>>17); ++k) list[i++] = k;
746         for (k = 4681 + (beg>>14); k <= 4681 + (end>>14); ++k) list[i++] = k;
747         return i;
748 }
749
750 ti_iter_t ti_iter_first()
751 {
752         ti_iter_t iter;
753         iter = calloc(1, sizeof(struct __ti_iter_t));
754         iter->from_first = 1;
755         return iter;
756 }
757
758 ti_iter_t ti_iter_query(const ti_index_t *idx, int tid, int beg, int end)
759 {
760         uint16_t *bins;
761         int i, n_bins, n_off;
762         pair64_t *off;
763         khint_t k;
764         khash_t(i) *index;
765         uint64_t min_off;
766         ti_iter_t iter = 0;
767
768         if (beg < 0) beg = 0;
769         if (end < beg) return 0;
770         // initialize the iterator
771         iter = calloc(1, sizeof(struct __ti_iter_t));
772         iter->idx = idx; iter->tid = tid; iter->beg = beg; iter->end = end; iter->i = -1;
773         // random access
774         bins = (uint16_t*)calloc(MAX_BIN, 2);
775         n_bins = reg2bins(beg, end, bins);
776         index = idx->index[tid];
777         if (idx->index2[tid].n > 0) {
778                 min_off = (beg>>TAD_LIDX_SHIFT >= idx->index2[tid].n)? idx->index2[tid].offset[idx->index2[tid].n-1]
779                         : idx->index2[tid].offset[beg>>TAD_LIDX_SHIFT];
780                 if (min_off == 0) { // improvement for index files built by tabix prior to 0.1.4
781                         int n = beg>>TAD_LIDX_SHIFT;
782                         if (n > idx->index2[tid].n) n = idx->index2[tid].n;
783                         for (i = n - 1; i >= 0; --i)
784                                 if (idx->index2[tid].offset[i] != 0) break;
785                         if (i >= 0) min_off = idx->index2[tid].offset[i];
786                 }
787         } else min_off = 0; // tabix 0.1.2 may produce such index files
788         for (i = n_off = 0; i < n_bins; ++i) {
789                 if ((k = kh_get(i, index, bins[i])) != kh_end(index))
790                         n_off += kh_value(index, k).n;
791         }
792         if (n_off == 0) {
793                 free(bins); return iter;
794         }
795         off = (pair64_t*)calloc(n_off, 16);
796         for (i = n_off = 0; i < n_bins; ++i) {
797                 if ((k = kh_get(i, index, bins[i])) != kh_end(index)) {
798                         int j;
799                         ti_binlist_t *p = &kh_value(index, k);
800                         for (j = 0; j < p->n; ++j)
801                                 if (p->list[j].v > min_off) off[n_off++] = p->list[j];
802                 }
803         }
804         free(bins);
805         {
806                 int l;
807                 ks_introsort(off, n_off, off);
808                 // resolve completely contained adjacent blocks
809                 for (i = 1, l = 0; i < n_off; ++i)
810                         if (off[l].v < off[i].v)
811                                 off[++l] = off[i];
812                 n_off = l + 1;
813                 // resolve overlaps between adjacent blocks; this may happen due to the merge in indexing
814                 for (i = 1; i < n_off; ++i)
815                         if (off[i-1].v >= off[i].u) off[i-1].v = off[i].u;
816                 { // merge adjacent blocks
817                         for (i = 1, l = 0; i < n_off; ++i) {
818                                 if (off[l].v>>16 == off[i].u>>16) off[l].v = off[i].v;
819                                 else off[++l] = off[i];
820                         }
821                         n_off = l + 1;
822                 }
823         }
824         iter->n_off = n_off; iter->off = off;
825         return iter;
826 }
827
828 const char *ti_iter_read(BGZF *fp, ti_iter_t iter, int *len)
829 {
830         if (iter->finished) return 0;
831         if (iter->from_first) {
832                 int ret;
833                 if ((ret = ti_readline(fp, &iter->str)) < 0) {
834                         iter->finished = 1;
835                         return 0;
836                 } else {
837                         if (len) *len = iter->str.l;
838                         return iter->str.s;
839                 }
840         }
841         if (iter->n_off == 0) return 0;
842         while (1) {
843                 int ret;
844                 if (iter->curr_off == 0 || iter->curr_off >= iter->off[iter->i].v) { // then jump to the next chunk
845                         if (iter->i == iter->n_off - 1) break; // no more chunks
846                         if (iter->i >= 0) assert(iter->curr_off == iter->off[iter->i].v); // otherwise bug
847                         if (iter->i < 0 || iter->off[iter->i].v != iter->off[iter->i+1].u) { // not adjacent chunks; then seek
848                                 bgzf_seek(fp, iter->off[iter->i+1].u, SEEK_SET);
849                                 iter->curr_off = bgzf_tell(fp);
850                         }
851                         ++iter->i;
852                 }
853                 if ((ret = ti_readline(fp, &iter->str)) >= 0) {
854                         ti_intv_t intv;
855                         iter->curr_off = bgzf_tell(fp);
856                         if (iter->str.s[0] == iter->idx->conf.meta_char) continue;
857                         get_intv((ti_index_t*)iter->idx, &iter->str, &intv);
858                         if (intv.tid != iter->tid || intv.beg >= iter->end) break; // no need to proceed
859                         else if (intv.end > iter->beg && iter->end > intv.beg) {
860                                 if (len) *len = iter->str.l;
861                                 return iter->str.s;
862                         }
863                 } else break; // end of file
864         }
865         iter->finished = 1;
866         return 0;
867 }
868
869 void ti_iter_destroy(ti_iter_t iter)
870 {
871         if (iter) {
872                 free(iter->str.s); free(iter->off);
873                 free(iter);
874         }
875 }
876
877 int ti_fetch(BGZF *fp, const ti_index_t *idx, int tid, int beg, int end, void *data, ti_fetch_f func)
878 {
879         ti_iter_t iter;
880         const char *s;
881         int len;
882         iter = ti_iter_query(idx, tid, beg, end);
883         while ((s = ti_iter_read(fp, iter, &len)) != 0)
884                 func(len, s, data);
885         ti_iter_destroy(iter);
886         return 0;
887 }
888
889 /*******************
890  * High-level APIs *
891  *******************/
892
893 tabix_t *ti_open(const char *fn, const char *fnidx)
894 {
895         tabix_t *t;
896         BGZF *fp;
897         if ((fp = bgzf_open(fn, "r")) == 0) return 0;
898         t = calloc(1, sizeof(tabix_t));
899         t->fn = strdup(fn);
900         if (fnidx) t->fnidx = strdup(fnidx);
901         t->fp = fp;
902         return t;
903 }
904
905 void ti_close(tabix_t *t)
906 {
907         if (t) {
908                 bgzf_close(t->fp);
909                 if (t->idx) ti_index_destroy(t->idx);
910                 free(t->fn); free(t->fnidx);
911                 free(t);
912         }
913 }
914
915 int ti_lazy_index_load(tabix_t *t)
916 {
917         if (t->idx == 0) { // load index
918                 if (t->fnidx) t->idx = ti_index_load_local(t->fnidx);
919                 else t->idx = ti_index_load(t->fn);
920                 if (t->idx == 0) return -1; // fail to load index
921         }
922         return 0;
923 }
924
925 ti_iter_t ti_queryi(tabix_t *t, int tid, int beg, int end)
926 {
927         if (tid < 0) return ti_iter_first();
928         if (ti_lazy_index_load(t) != 0) return 0;
929         return ti_iter_query(t->idx, tid, beg, end);    
930 }
931
932 ti_iter_t ti_query(tabix_t *t, const char *name, int beg, int end)
933 {
934         int tid;
935         if (name == 0) return ti_iter_first();
936         // then need to load the index
937         if (ti_lazy_index_load(t) != 0) return 0;
938         if ((tid = ti_get_tid(t->idx, name)) < 0) return 0;
939         return ti_iter_query(t->idx, tid, beg, end);
940 }
941
942 const char *ti_read(tabix_t *t, ti_iter_t iter, int *len)
943 {
944         return ti_iter_read(t->fp, iter, len);
945 }