# New upstream release; sent ITP.
[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(offt, 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 == 4) {
196                                                 if (b < i) intv->end = intv->beg + (i - b);
197                                         } else if (id == 8) { // look for "END="
198                                                 int c = str->s[i];
199                                                 str->s[i] = 0;
200                                                 s = strstr(str->s + b, "END=");
201                                                 if (s == str->s + b) s += 4;
202                                                 else if (s) {
203                                                         s = strstr(str->s + b, ";END=");
204                                                         if (s) s += 5;
205                                                 }
206                                                 if (s) intv->end = strtol(s, &s, 0);
207                                                 str->s[i] = c;
208                                         }
209                                 }
210                         }
211                         b = i + 1;
212                         ++id;
213                 }
214         }
215         if (intv->tid < 0 || intv->beg < 0 || intv->end < 0) return -1;
216         intv->bin = ti_reg2bin(intv->beg, intv->end);
217         return 0;
218 }
219
220 /************
221  * indexing *
222  ************/
223
224 // requirement: len <= LEN_MASK
225 static inline void insert_offset(khash_t(i) *h, int bin, uint64_t beg, uint64_t end)
226 {
227         khint_t k;
228         ti_binlist_t *l;
229         int ret;
230         k = kh_put(i, h, bin, &ret);
231         l = &kh_value(h, k);
232         if (ret) { // not present
233                 l->m = 1; l->n = 0;
234                 l->list = (pair64_t*)calloc(l->m, 16);
235         }
236         if (l->n == l->m) {
237                 l->m <<= 1;
238                 l->list = (pair64_t*)realloc(l->list, l->m * 16);
239         }
240         l->list[l->n].u = beg; l->list[l->n++].v = end;
241 }
242
243 static inline uint64_t insert_offset2(ti_lidx_t *index2, int _beg, int _end, uint64_t offset)
244 {
245         int i, beg, end;
246         beg = _beg >> TAD_LIDX_SHIFT;
247         end = (_end - 1) >> TAD_LIDX_SHIFT;
248         if (index2->m < end + 1) {
249                 int old_m = index2->m;
250                 index2->m = end + 1;
251                 kroundup32(index2->m);
252                 index2->offset = (uint64_t*)realloc(index2->offset, index2->m * 8);
253                 memset(index2->offset + old_m, 0, 8 * (index2->m - old_m));
254         }
255         if (beg == end) {
256                 if (index2->offset[beg] == 0) index2->offset[beg] = offset;
257         } else {
258                 for (i = beg; i <= end; ++i)
259                         if (index2->offset[i] == 0) index2->offset[i] = offset;
260         }
261         if (index2->n < end + 1) index2->n = end + 1;
262         return (uint64_t)beg<<32 | end;
263 }
264
265 static void merge_chunks(ti_index_t *idx)
266 {
267         khash_t(i) *index;
268         int i, l, m;
269         khint_t k;
270         for (i = 0; i < idx->n; ++i) {
271                 index = idx->index[i];
272                 for (k = kh_begin(index); k != kh_end(index); ++k) {
273                         ti_binlist_t *p;
274                         if (!kh_exist(index, k)) continue;
275                         p = &kh_value(index, k);
276                         m = 0;
277                         for (l = 1; l < p->n; ++l) {
278                                 if (p->list[m].v>>16 == p->list[l].u>>16) p->list[m].v = p->list[l].v;
279                                 else p->list[++m] = p->list[l];
280                         } // ~for(l)
281                         p->n = m + 1;
282                 } // ~for(k)
283         } // ~for(i)
284 }
285
286 static void fill_missing(ti_index_t *idx)
287 {
288         int i, j;
289         for (i = 0; i < idx->n; ++i) {
290                 ti_lidx_t *idx2 = &idx->index2[i];
291                 for (j = 1; j < idx2->n; ++j)
292                         if (idx2->offset[j] == 0)
293                                 idx2->offset[j] = idx2->offset[j-1];
294         }
295 }
296
297 ti_index_t *ti_index_core(BGZF *fp, const ti_conf_t *conf)
298 {
299         int ret;
300         ti_index_t *idx;
301         uint32_t last_bin, save_bin;
302         int32_t last_coor, last_tid, save_tid;
303         uint64_t save_off, last_off, lineno = 0, offset0 = (uint64_t)-1, tmp;
304         kstring_t *str;
305
306         str = calloc(1, sizeof(kstring_t));
307
308         idx = (ti_index_t*)calloc(1, sizeof(ti_index_t));
309         idx->conf = *conf;
310         idx->n = idx->max = 0;
311         idx->tname = kh_init(s);
312         idx->index = 0;
313         idx->index2 = 0;
314
315         save_bin = save_tid = last_tid = last_bin = 0xffffffffu;
316         save_off = last_off = bgzf_tell(fp); last_coor = 0xffffffffu;
317         while ((ret = ti_readline(fp, str)) >= 0) {
318                 ti_intv_t intv;
319                 ++lineno;
320                 if (lineno <= idx->conf.line_skip || str->s[0] == idx->conf.meta_char) {
321                         last_off = bgzf_tell(fp);
322                         continue;
323                 }
324                 get_intv(idx, str, &intv);
325                 if (last_tid != intv.tid) { // change of chromosomes
326             if (last_tid>intv.tid )
327             {
328                 fprintf(stderr,"[ti_index_core] the chromosome blocks not continuous at line %llu, is the file sorted?\n",(unsigned long long)lineno);
329                 exit(1);
330             }
331                         last_tid = intv.tid;
332                         last_bin = 0xffffffffu;
333                 } else if (last_coor > intv.beg) {
334                         fprintf(stderr, "[ti_index_core] the file out of order at line %llu\n", (unsigned long long)lineno);
335                         exit(1);
336                 }
337                 tmp = insert_offset2(&idx->index2[intv.tid], intv.beg, intv.end, last_off);
338                 if (last_off == 0) offset0 = tmp;
339                 if (intv.bin != last_bin) { // then possibly write the binning index
340                         if (save_bin != 0xffffffffu) // save_bin==0xffffffffu only happens to the first record
341                                 insert_offset(idx->index[save_tid], save_bin, save_off, last_off);
342                         save_off = last_off;
343                         save_bin = last_bin = intv.bin;
344                         save_tid = intv.tid;
345                         if (save_tid < 0) break;
346                 }
347                 if (bgzf_tell(fp) <= last_off) {
348                         fprintf(stderr, "[ti_index_core] bug in BGZF: %llx < %llx\n",
349                                         (unsigned long long)bgzf_tell(fp), (unsigned long long)last_off);
350                         exit(1);
351                 }
352                 last_off = bgzf_tell(fp);
353                 last_coor = intv.beg;
354         }
355         if (save_tid >= 0) insert_offset(idx->index[save_tid], save_bin, save_off, bgzf_tell(fp));
356         merge_chunks(idx);
357         fill_missing(idx);
358         if (offset0 != (uint64_t)-1 && idx->n && idx->index2[0].offset) {
359                 int i, beg = offset0>>32, end = offset0&0xffffffffu;
360                 for (i = beg; i <= end; ++i) idx->index2[0].offset[i] = 0;
361         }
362
363         free(str->s); free(str);
364         return idx;
365 }
366
367 void ti_index_destroy(ti_index_t *idx)
368 {
369         khint_t k;
370         int i;
371         if (idx == 0) return;
372         // destroy the name hash table
373         for (k = kh_begin(idx->tname); k != kh_end(idx->tname); ++k) {
374                 if (kh_exist(idx->tname, k))
375                         free((char*)kh_key(idx->tname, k));
376         }
377         kh_destroy(s, idx->tname);
378         // destroy the binning index
379         for (i = 0; i < idx->n; ++i) {
380                 khash_t(i) *index = idx->index[i];
381                 ti_lidx_t *index2 = idx->index2 + i;
382                 for (k = kh_begin(index); k != kh_end(index); ++k) {
383                         if (kh_exist(index, k))
384                                 free(kh_value(index, k).list);
385                 }
386                 kh_destroy(i, index);
387                 free(index2->offset);
388         }
389         free(idx->index);
390         // destroy the linear index
391         free(idx->index2);
392         free(idx);
393 }
394
395 /******************
396  * index file I/O *
397  ******************/
398
399 void ti_index_save(const ti_index_t *idx, BGZF *fp)
400 {
401         int32_t i, size, ti_is_be;
402         khint_t k;
403         ti_is_be = bam_is_big_endian();
404         bgzf_write(fp, "TBI\1", 4);
405         if (ti_is_be) {
406                 uint32_t x = idx->n;
407                 bgzf_write(fp, bam_swap_endian_4p(&x), 4);
408         } else bgzf_write(fp, &idx->n, 4);
409         assert(sizeof(ti_conf_t) == 24);
410         if (ti_is_be) { // write ti_conf_t;
411                 uint32_t x[6];
412                 memcpy(x, &idx->conf, 24);
413                 for (i = 0; i < 6; ++i) bgzf_write(fp, bam_swap_endian_4p(&x[i]), 4);
414         } else bgzf_write(fp, &idx->conf, sizeof(ti_conf_t));
415         { // write target names
416                 char **name;
417                 int32_t l = 0;
418                 name = calloc(kh_size(idx->tname), sizeof(void*));
419                 for (k = kh_begin(idx->tname); k != kh_end(idx->tname); ++k)
420                         if (kh_exist(idx->tname, k))
421                                 name[kh_value(idx->tname, k)] = (char*)kh_key(idx->tname, k);
422                 for (i = 0; i < kh_size(idx->tname); ++i)
423                         l += strlen(name[i]) + 1;
424                 if (ti_is_be) bgzf_write(fp, bam_swap_endian_4p(&l), 4);
425                 else bgzf_write(fp, &l, 4);
426                 for (i = 0; i < kh_size(idx->tname); ++i)
427                         bgzf_write(fp, name[i], strlen(name[i]) + 1);
428                 free(name);
429         }
430         for (i = 0; i < idx->n; ++i) {
431                 khash_t(i) *index = idx->index[i];
432                 ti_lidx_t *index2 = idx->index2 + i;
433                 // write binning index
434                 size = kh_size(index);
435                 if (ti_is_be) { // big endian
436                         uint32_t x = size;
437                         bgzf_write(fp, bam_swap_endian_4p(&x), 4);
438                 } else bgzf_write(fp, &size, 4);
439                 for (k = kh_begin(index); k != kh_end(index); ++k) {
440                         if (kh_exist(index, k)) {
441                                 ti_binlist_t *p = &kh_value(index, k);
442                                 if (ti_is_be) { // big endian
443                                         uint32_t x;
444                                         x = kh_key(index, k); bgzf_write(fp, bam_swap_endian_4p(&x), 4);
445                                         x = p->n; bgzf_write(fp, bam_swap_endian_4p(&x), 4);
446                                         for (x = 0; (int)x < p->n; ++x) {
447                                                 bam_swap_endian_8p(&p->list[x].u);
448                                                 bam_swap_endian_8p(&p->list[x].v);
449                                         }
450                                         bgzf_write(fp, p->list, 16 * p->n);
451                                         for (x = 0; (int)x < p->n; ++x) {
452                                                 bam_swap_endian_8p(&p->list[x].u);
453                                                 bam_swap_endian_8p(&p->list[x].v);
454                                         }
455                                 } else {
456                                         bgzf_write(fp, &kh_key(index, k), 4);
457                                         bgzf_write(fp, &p->n, 4);
458                                         bgzf_write(fp, p->list, 16 * p->n);
459                                 }
460                         }
461                 }
462                 // write linear index (index2)
463                 if (ti_is_be) {
464                         int x = index2->n;
465                         bgzf_write(fp, bam_swap_endian_4p(&x), 4);
466                 } else bgzf_write(fp, &index2->n, 4);
467                 if (ti_is_be) { // big endian
468                         int x;
469                         for (x = 0; (int)x < index2->n; ++x)
470                                 bam_swap_endian_8p(&index2->offset[x]);
471                         bgzf_write(fp, index2->offset, 8 * index2->n);
472                         for (x = 0; (int)x < index2->n; ++x)
473                                 bam_swap_endian_8p(&index2->offset[x]);
474                 } else bgzf_write(fp, index2->offset, 8 * index2->n);
475         }
476 }
477
478 static ti_index_t *ti_index_load_core(BGZF *fp)
479 {
480         int i, ti_is_be;
481         char magic[4];
482         ti_index_t *idx;
483         ti_is_be = bam_is_big_endian();
484         if (fp == 0) {
485                 fprintf(stderr, "[ti_index_load_core] fail to load index.\n");
486                 return 0;
487         }
488         bgzf_read(fp, magic, 4);
489         if (strncmp(magic, "TBI\1", 4)) {
490                 fprintf(stderr, "[ti_index_load] wrong magic number.\n");
491                 return 0;
492         }
493         idx = (ti_index_t*)calloc(1, sizeof(ti_index_t));       
494         bgzf_read(fp, &idx->n, 4);
495         if (ti_is_be) bam_swap_endian_4p(&idx->n);
496         idx->tname = kh_init(s);
497         idx->index = (khash_t(i)**)calloc(idx->n, sizeof(void*));
498         idx->index2 = (ti_lidx_t*)calloc(idx->n, sizeof(ti_lidx_t));
499         // read idx->conf
500         bgzf_read(fp, &idx->conf, sizeof(ti_conf_t));
501         if (ti_is_be) {
502                 bam_swap_endian_4p(&idx->conf.preset);
503                 bam_swap_endian_4p(&idx->conf.sc);
504                 bam_swap_endian_4p(&idx->conf.bc);
505                 bam_swap_endian_4p(&idx->conf.ec);
506                 bam_swap_endian_4p(&idx->conf.meta_char);
507                 bam_swap_endian_4p(&idx->conf.line_skip);
508         }
509         { // read target names
510                 int j, ret;
511                 kstring_t *str;
512                 int32_t l;
513                 uint8_t *buf;
514                 bgzf_read(fp, &l, 4);
515                 if (ti_is_be) bam_swap_endian_4p(&l);
516                 buf = calloc(l, 1);
517                 bgzf_read(fp, buf, l);
518                 str = calloc(1, sizeof(kstring_t));
519                 for (i = j = 0; i < l; ++i) {
520                         if (buf[i] == 0) {
521                                 khint_t k = kh_put(s, idx->tname, strdup(str->s), &ret);
522                                 kh_value(idx->tname, k) = j++;
523                                 str->l = 0;
524                         } else kputc(buf[i], str);
525                 }
526                 free(str->s); free(str); free(buf);
527         }
528         for (i = 0; i < idx->n; ++i) {
529                 khash_t(i) *index;
530                 ti_lidx_t *index2 = idx->index2 + i;
531                 uint32_t key, size;
532                 khint_t k;
533                 int j, ret;
534                 ti_binlist_t *p;
535                 index = idx->index[i] = kh_init(i);
536                 // load binning index
537                 bgzf_read(fp, &size, 4);
538                 if (ti_is_be) bam_swap_endian_4p(&size);
539                 for (j = 0; j < (int)size; ++j) {
540                         bgzf_read(fp, &key, 4);
541                         if (ti_is_be) bam_swap_endian_4p(&key);
542                         k = kh_put(i, index, key, &ret);
543                         p = &kh_value(index, k);
544                         bgzf_read(fp, &p->n, 4);
545                         if (ti_is_be) bam_swap_endian_4p(&p->n);
546                         p->m = p->n;
547                         p->list = (pair64_t*)malloc(p->m * 16);
548                         bgzf_read(fp, p->list, 16 * p->n);
549                         if (ti_is_be) {
550                                 int x;
551                                 for (x = 0; x < p->n; ++x) {
552                                         bam_swap_endian_8p(&p->list[x].u);
553                                         bam_swap_endian_8p(&p->list[x].v);
554                                 }
555                         }
556                 }
557                 // load linear index
558                 bgzf_read(fp, &index2->n, 4);
559                 if (ti_is_be) bam_swap_endian_4p(&index2->n);
560                 index2->m = index2->n;
561                 index2->offset = (uint64_t*)calloc(index2->m, 8);
562                 bgzf_read(fp, index2->offset, index2->n * 8);
563                 if (ti_is_be)
564                         for (j = 0; j < index2->n; ++j) bam_swap_endian_8p(&index2->offset[j]);
565         }
566         return idx;
567 }
568
569 ti_index_t *ti_index_load_local(const char *fnidx)
570 {
571         BGZF *fp;
572         fp = bgzf_open(fnidx, "r");
573         if (fp) {
574                 ti_index_t *idx = ti_index_load_core(fp);
575                 bgzf_close(fp);
576                 return idx;
577         } else return 0;
578 }
579
580 #ifdef _USE_KNETFILE
581 static void download_from_remote(const char *url)
582 {
583         const int buf_size = 1 * 1024 * 1024;
584         char *fn;
585         FILE *fp;
586         uint8_t *buf;
587         knetFile *fp_remote;
588         int l;
589         if (strstr(url, "ftp://") != url && strstr(url, "http://") != url) return;
590         l = strlen(url);
591         for (fn = (char*)url + l - 1; fn >= url; --fn)
592                 if (*fn == '/') break;
593         ++fn; // fn now points to the file name
594         fp_remote = knet_open(url, "r");
595         if (fp_remote == 0) {
596                 fprintf(stderr, "[download_from_remote] fail to open remote file.\n");
597                 return;
598         }
599         if ((fp = fopen(fn, "w")) == 0) {
600                 fprintf(stderr, "[download_from_remote] fail to create file in the working directory.\n");
601                 knet_close(fp_remote);
602                 return;
603         }
604         buf = (uint8_t*)calloc(buf_size, 1);
605         while ((l = knet_read(fp_remote, buf, buf_size)) != 0)
606                 fwrite(buf, 1, l, fp);
607         free(buf);
608         fclose(fp);
609         knet_close(fp_remote);
610 }
611 #else
612 static void download_from_remote(const char *url)
613 {
614         return;
615 }
616 #endif
617
618 static char *get_local_version(const char *fn)
619 {
620     struct stat sbuf;
621         char *fnidx = (char*)calloc(strlen(fn) + 5, 1);
622         strcat(strcpy(fnidx, fn), ".tbi");
623         if ((strstr(fnidx, "ftp://") == fnidx || strstr(fnidx, "http://") == fnidx)) {
624                 char *p, *url;
625                 int l = strlen(fnidx);
626                 for (p = fnidx + l - 1; p >= fnidx; --p)
627                         if (*p == '/') break;
628                 url = fnidx; fnidx = strdup(p + 1);
629                 if (stat(fnidx, &sbuf) == 0) {
630                         free(url);
631                         return fnidx;
632                 }
633                 fprintf(stderr, "[%s] downloading the index file...\n", __func__);
634                 download_from_remote(url);
635                 free(url);
636         }
637     if (stat(fnidx, &sbuf) == 0) return fnidx;
638         free(fnidx); return 0;
639 }
640
641 const char **ti_seqname(const ti_index_t *idx, int *n)
642 {
643         const char **names;
644         khint_t k;
645         *n = idx->n;
646         names = calloc(idx->n, sizeof(void*));
647         for (k = kh_begin(idx->tname); k < kh_end(idx->tname); ++k)
648                 if (kh_exist(idx->tname, k))
649                         names[kh_val(idx->tname, k)] = kh_key(idx->tname, k);
650         return names;
651 }
652
653 ti_index_t *ti_index_load(const char *fn)
654 {
655         ti_index_t *idx;
656     char *fname = get_local_version(fn);
657         if (fname == 0) return 0;
658         idx = ti_index_load_local(fname);
659         if (idx == 0) fprintf(stderr, "[ti_index_load] fail to load the index: %s\n", fname);
660     free(fname);
661         return idx;
662 }
663
664 int ti_index_build2(const char *fn, const ti_conf_t *conf, const char *_fnidx)
665 {
666         char *fnidx;
667         BGZF *fp, *fpidx;
668         ti_index_t *idx;
669         if ((fp = bgzf_open(fn, "r")) == 0) {
670                 fprintf(stderr, "[ti_index_build2] fail to open the file: %s\n", fn);
671                 return -1;
672         }
673         idx = ti_index_core(fp, conf);
674         bgzf_close(fp);
675         if (_fnidx == 0) {
676                 fnidx = (char*)calloc(strlen(fn) + 5, 1);
677                 strcpy(fnidx, fn); strcat(fnidx, ".tbi");
678         } else fnidx = strdup(_fnidx);
679         fpidx = bgzf_open(fnidx, "w");
680         if (fpidx == 0) {
681                 fprintf(stderr, "[ti_index_build2] fail to create the index file.\n");
682                 free(fnidx);
683                 return -1;
684         }
685         ti_index_save(idx, fpidx);
686         ti_index_destroy(idx);
687         bgzf_close(fpidx);
688         free(fnidx);
689         return 0;
690 }
691
692 int ti_index_build(const char *fn, const ti_conf_t *conf)
693 {
694         return ti_index_build2(fn, conf, 0);
695 }
696
697 /********************************************
698  * parse a region in the format chr:beg-end *
699  ********************************************/
700
701 int ti_get_tid(const ti_index_t *idx, const char *name)
702 {
703         khiter_t iter;
704         const khash_t(s) *h = idx->tname;
705         iter = kh_get(s, h, name); /* get the tid */
706         if (iter == kh_end(h)) return -1;
707         return kh_value(h, iter);
708 }
709
710 int ti_parse_region(const ti_index_t *idx, const char *str, int *tid, int *begin, int *end)
711 {
712         char *s, *p;
713         int i, l, k;
714         l = strlen(str);
715         p = s = (char*)malloc(l+1);
716         /* squeeze out "," */
717         for (i = k = 0; i != l; ++i)
718                 if (str[i] != ',' && !isspace(str[i])) s[k++] = str[i];
719         s[k] = 0;
720         for (i = 0; i != k; ++i) if (s[i] == ':') break;
721         s[i] = 0;
722         if ((*tid = ti_get_tid(idx, s)) < 0) {
723                 free(s);
724                 return -1;
725         }
726         if (i == k) { /* dump the whole sequence */
727                 *begin = 0; *end = 1<<29; free(s);
728                 return 0;
729         }
730         for (p = s + i + 1; i != k; ++i) if (s[i] == '-') break;
731         *begin = atoi(p);
732         if (i < k) {
733                 p = s + i + 1;
734                 *end = atoi(p);
735         } else *end = 1<<29;
736         if (*begin > 0) --*begin;
737         free(s);
738         if (*begin > *end) return -1;
739         return 0;
740 }
741
742 /*******************************
743  * retrieve a specified region *
744  *******************************/
745
746 #define MAX_BIN 37450 // =(8^6-1)/7+1
747
748 static inline int reg2bins(uint32_t beg, uint32_t end, uint16_t list[MAX_BIN])
749 {
750         int i = 0, k;
751         if (beg >= end) return 0;
752         if (end >= 1u<<29) end = 1u<<29;
753         --end;
754         list[i++] = 0;
755         for (k =    1 + (beg>>26); k <=    1 + (end>>26); ++k) list[i++] = k;
756         for (k =    9 + (beg>>23); k <=    9 + (end>>23); ++k) list[i++] = k;
757         for (k =   73 + (beg>>20); k <=   73 + (end>>20); ++k) list[i++] = k;
758         for (k =  585 + (beg>>17); k <=  585 + (end>>17); ++k) list[i++] = k;
759         for (k = 4681 + (beg>>14); k <= 4681 + (end>>14); ++k) list[i++] = k;
760         return i;
761 }
762
763 ti_iter_t ti_iter_first()
764 {
765         ti_iter_t iter;
766         iter = calloc(1, sizeof(struct __ti_iter_t));
767         iter->from_first = 1;
768         return iter;
769 }
770
771 ti_iter_t ti_iter_query(const ti_index_t *idx, int tid, int beg, int end)
772 {
773         uint16_t *bins;
774         int i, n_bins, n_off;
775         pair64_t *off;
776         khint_t k;
777         khash_t(i) *index;
778         uint64_t min_off;
779         ti_iter_t iter = 0;
780
781         if (beg < 0) beg = 0;
782         if (end < beg) return 0;
783         // initialize the iterator
784         iter = calloc(1, sizeof(struct __ti_iter_t));
785         iter->idx = idx; iter->tid = tid; iter->beg = beg; iter->end = end; iter->i = -1;
786         // random access
787         bins = (uint16_t*)calloc(MAX_BIN, 2);
788         n_bins = reg2bins(beg, end, bins);
789         index = idx->index[tid];
790         if (idx->index2[tid].n > 0) {
791                 min_off = (beg>>TAD_LIDX_SHIFT >= idx->index2[tid].n)? idx->index2[tid].offset[idx->index2[tid].n-1]
792                         : idx->index2[tid].offset[beg>>TAD_LIDX_SHIFT];
793                 if (min_off == 0) { // improvement for index files built by tabix prior to 0.1.4
794                         int n = beg>>TAD_LIDX_SHIFT;
795                         if (n > idx->index2[tid].n) n = idx->index2[tid].n;
796                         for (i = n - 1; i >= 0; --i)
797                                 if (idx->index2[tid].offset[i] != 0) break;
798                         if (i >= 0) min_off = idx->index2[tid].offset[i];
799                 }
800         } else min_off = 0; // tabix 0.1.2 may produce such index files
801         for (i = n_off = 0; i < n_bins; ++i) {
802                 if ((k = kh_get(i, index, bins[i])) != kh_end(index))
803                         n_off += kh_value(index, k).n;
804         }
805         if (n_off == 0) {
806                 free(bins); return iter;
807         }
808         off = (pair64_t*)calloc(n_off, 16);
809         for (i = n_off = 0; i < n_bins; ++i) {
810                 if ((k = kh_get(i, index, bins[i])) != kh_end(index)) {
811                         int j;
812                         ti_binlist_t *p = &kh_value(index, k);
813                         for (j = 0; j < p->n; ++j)
814                                 if (p->list[j].v > min_off) off[n_off++] = p->list[j];
815                 }
816         }
817         if (n_off == 0) {
818                 free(bins); free(off); return iter;
819         }
820         free(bins);
821         {
822                 int l;
823                 ks_introsort(offt, n_off, off);
824                 // resolve completely contained adjacent blocks
825                 for (i = 1, l = 0; i < n_off; ++i)
826                         if (off[l].v < off[i].v)
827                                 off[++l] = off[i];
828                 n_off = l + 1;
829                 // resolve overlaps between adjacent blocks; this may happen due to the merge in indexing
830                 for (i = 1; i < n_off; ++i)
831                         if (off[i-1].v >= off[i].u) off[i-1].v = off[i].u;
832                 { // merge adjacent blocks
833                         for (i = 1, l = 0; i < n_off; ++i) {
834                                 if (off[l].v>>16 == off[i].u>>16) off[l].v = off[i].v;
835                                 else off[++l] = off[i];
836                         }
837                         n_off = l + 1;
838                 }
839         }
840         iter->n_off = n_off; iter->off = off;
841         return iter;
842 }
843
844 const char *ti_iter_read(BGZF *fp, ti_iter_t iter, int *len)
845 {
846         if (iter->finished) return 0;
847         if (iter->from_first) {
848                 int ret;
849                 if ((ret = ti_readline(fp, &iter->str)) < 0) {
850                         iter->finished = 1;
851                         return 0;
852                 } else {
853                         if (len) *len = iter->str.l;
854                         return iter->str.s;
855                 }
856         }
857         if (iter->n_off == 0) return 0;
858         while (1) {
859                 int ret;
860                 if (iter->curr_off == 0 || iter->curr_off >= iter->off[iter->i].v) { // then jump to the next chunk
861                         if (iter->i == iter->n_off - 1) break; // no more chunks
862                         if (iter->i >= 0) assert(iter->curr_off == iter->off[iter->i].v); // otherwise bug
863                         if (iter->i < 0 || iter->off[iter->i].v != iter->off[iter->i+1].u) { // not adjacent chunks; then seek
864                                 bgzf_seek(fp, iter->off[iter->i+1].u, SEEK_SET);
865                                 iter->curr_off = bgzf_tell(fp);
866                         }
867                         ++iter->i;
868                 }
869                 if ((ret = ti_readline(fp, &iter->str)) >= 0) {
870                         ti_intv_t intv;
871                         iter->curr_off = bgzf_tell(fp);
872                         if (iter->str.s[0] == iter->idx->conf.meta_char) continue;
873                         get_intv((ti_index_t*)iter->idx, &iter->str, &intv);
874                         if (intv.tid != iter->tid || intv.beg >= iter->end) break; // no need to proceed
875                         else if (intv.end > iter->beg && iter->end > intv.beg) {
876                                 if (len) *len = iter->str.l;
877                                 return iter->str.s;
878                         }
879                 } else break; // end of file
880         }
881         iter->finished = 1;
882         return 0;
883 }
884
885 void ti_iter_destroy(ti_iter_t iter)
886 {
887         if (iter) {
888                 free(iter->str.s); free(iter->off);
889                 free(iter);
890         }
891 }
892
893 int ti_fetch(BGZF *fp, const ti_index_t *idx, int tid, int beg, int end, void *data, ti_fetch_f func)
894 {
895         ti_iter_t iter;
896         const char *s;
897         int len;
898         iter = ti_iter_query(idx, tid, beg, end);
899         while ((s = ti_iter_read(fp, iter, &len)) != 0)
900                 func(len, s, data);
901         ti_iter_destroy(iter);
902         return 0;
903 }
904
905 /*******************
906  * High-level APIs *
907  *******************/
908
909 tabix_t *ti_open(const char *fn, const char *fnidx)
910 {
911         tabix_t *t;
912         BGZF *fp;
913         if ((fp = bgzf_open(fn, "r")) == 0) return 0;
914         t = calloc(1, sizeof(tabix_t));
915         t->fn = strdup(fn);
916         if (fnidx) t->fnidx = strdup(fnidx);
917         t->fp = fp;
918         return t;
919 }
920
921 void ti_close(tabix_t *t)
922 {
923         if (t) {
924                 bgzf_close(t->fp);
925                 if (t->idx) ti_index_destroy(t->idx);
926                 free(t->fn); free(t->fnidx);
927                 free(t);
928         }
929 }
930
931 int ti_lazy_index_load(tabix_t *t)
932 {
933         if (t->idx == 0) { // load index
934                 if (t->fnidx) t->idx = ti_index_load_local(t->fnidx);
935                 else t->idx = ti_index_load(t->fn);
936                 if (t->idx == 0) return -1; // fail to load index
937         }
938         return 0;
939 }
940
941 ti_iter_t ti_queryi(tabix_t *t, int tid, int beg, int end)
942 {
943         if (tid < 0) return ti_iter_first();
944         if (ti_lazy_index_load(t) != 0) return 0;
945         return ti_iter_query(t->idx, tid, beg, end);    
946 }
947
948 ti_iter_t ti_querys(tabix_t *t, const char *reg)
949 {
950         int tid, beg, end;
951         if (reg == 0) return ti_iter_first();
952         if (ti_lazy_index_load(t) != 0) return 0;
953         if (ti_parse_region(t->idx, reg, &tid, &beg, &end) < 0) return 0;
954         return ti_iter_query(t->idx, tid, beg, end);
955 }
956
957 ti_iter_t ti_query(tabix_t *t, const char *name, int beg, int end)
958 {
959         int tid;
960         if (name == 0) return ti_iter_first();
961         // then need to load the index
962         if (ti_lazy_index_load(t) != 0) return 0;
963         if ((tid = ti_get_tid(t->idx, name)) < 0) return 0;
964         return ti_iter_query(t->idx, tid, beg, end);
965 }
966
967 const char *ti_read(tabix_t *t, ti_iter_t iter, int *len)
968 {
969         return ti_iter_read(t->fp, iter, len);
970 }