7 #include "bam_endian.h"
13 #define TAD_MIN_CHUNK_GAP 32768
14 // 1<<14 is the size of minimum bin.
15 #define TAD_LIDX_SHIFT 14
21 #define pair64_lt(a,b) ((a).u < (b).u)
22 KSORT_INIT(offt, pair64_t, pair64_lt)
34 KHASH_MAP_INIT_INT(i, ti_binlist_t)
35 KHASH_MAP_INIT_STR(s, int)
46 int from_first; // read from the first record; no random access
47 int tid, beg, end, n_off, i, finished;
50 const ti_index_t *idx;
55 int tid, beg, end, bin;
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 };
69 int ti_readline(BGZF *fp, kstring_t *str)
73 while ((c = bgzf_getc(fp)) >= 0 && c != '\n') {
75 if (c != '\r') kputc(c, str);
77 if (c < 0 && l == 0) return -1; // end of file
82 /* Below is a faster implementation largely equivalent to the one
83 * commented out above. */
84 int ti_readline(BGZF *fp, kstring_t *str)
87 unsigned char *buf = (unsigned char*)fp->uncompressed_block;
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; }
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;
100 str->s = (char*)realloc(str->s, str->m);
102 memcpy(str->s + str->l, buf + fp->block_offset, l);
104 fp->block_offset += l + 1;
105 if (fp->block_offset >= fp->block_length) {
107 fp->block_address = knet_tell(fp->x.fpr);
109 fp->block_address = ftello(fp->file);
111 fp->block_offset = 0;
112 fp->block_length = 0;
114 } while (state == 0);
115 if (str->l == 0 && state < 0) return state;
120 /*************************************
121 * get the interval from a data line *
122 *************************************/
124 static inline int ti_reg2bin(uint32_t beg, uint32_t 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);
135 static int get_tid(ti_index_t *idx, const char *ss)
139 k = kh_get(s, idx->tname, ss);
140 if (k == kh_end(idx->tname)) { // a new target sequence
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));
148 memset(&idx->index2[idx->n], 0, sizeof(ti_lidx_t));
149 idx->index[idx->n++] = kh_init(i);
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);
159 static int get_intv(ti_index_t *idx, kstring_t *str, ti_intv_t *intv)
161 int i, b = 0, id = 1;
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) {
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;
175 if (intv->beg < 0) intv->beg = 0;
176 if (intv->end < 1) intv->end = 1;
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
184 for (s = str->s + b; s < str->s + i;) {
185 long x = strtol(s, &t, 10);
187 if (op == 'M' || op == 'D' || op == 'N') l += x;
191 intv->end = intv->beg + l;
193 } else if ((idx->conf.preset&0xffff) == TI_PRESET_VCF) {
194 // FIXME: the following is NOT tested and is likely to be buggy
196 if (b < i) intv->end = intv->beg + (i - b);
197 } else if (id == 8) { // look for "END="
200 s = strstr(str->s + b, "END=");
201 if (s == str->s + b) s += 4;
203 s = strstr(str->s + b, ";END=");
206 if (s) intv->end = strtol(s, &s, 0);
215 if (intv->tid < 0 || intv->beg < 0 || intv->end < 0) return -1;
216 intv->bin = ti_reg2bin(intv->beg, intv->end);
224 // requirement: len <= LEN_MASK
225 static inline void insert_offset(khash_t(i) *h, int bin, uint64_t beg, uint64_t end)
230 k = kh_put(i, h, bin, &ret);
232 if (ret) { // not present
234 l->list = (pair64_t*)calloc(l->m, 16);
238 l->list = (pair64_t*)realloc(l->list, l->m * 16);
240 l->list[l->n].u = beg; l->list[l->n++].v = end;
243 static inline uint64_t insert_offset2(ti_lidx_t *index2, int _beg, int _end, uint64_t offset)
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;
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));
256 if (index2->offset[beg] == 0) index2->offset[beg] = offset;
258 for (i = beg; i <= end; ++i)
259 if (index2->offset[i] == 0) index2->offset[i] = offset;
261 if (index2->n < end + 1) index2->n = end + 1;
262 return (uint64_t)beg<<32 | end;
265 static void merge_chunks(ti_index_t *idx)
270 for (i = 0; i < idx->n; ++i) {
271 index = idx->index[i];
272 for (k = kh_begin(index); k != kh_end(index); ++k) {
274 if (!kh_exist(index, k)) continue;
275 p = &kh_value(index, k);
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];
286 static void fill_missing(ti_index_t *idx)
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];
297 ti_index_t *ti_index_core(BGZF *fp, const ti_conf_t *conf)
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;
306 str = calloc(1, sizeof(kstring_t));
308 idx = (ti_index_t*)calloc(1, sizeof(ti_index_t));
310 idx->n = idx->max = 0;
311 idx->tname = kh_init(s);
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) {
320 if (lineno <= idx->conf.line_skip || str->s[0] == idx->conf.meta_char) {
321 last_off = bgzf_tell(fp);
324 get_intv(idx, str, &intv);
325 if (last_tid != intv.tid) { // change of chromosomes
326 if (last_tid>intv.tid )
328 fprintf(stderr,"[ti_index_core] the chromosome blocks not continuous at line %llu, is the file sorted?\n",(unsigned long long)lineno);
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);
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);
343 save_bin = last_bin = intv.bin;
345 if (save_tid < 0) break;
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);
352 last_off = bgzf_tell(fp);
353 last_coor = intv.beg;
355 if (save_tid >= 0) insert_offset(idx->index[save_tid], save_bin, save_off, bgzf_tell(fp));
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;
363 free(str->s); free(str);
367 void ti_index_destroy(ti_index_t *idx)
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));
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);
386 kh_destroy(i, index);
387 free(index2->offset);
390 // destroy the linear index
399 void ti_index_save(const ti_index_t *idx, BGZF *fp)
401 int32_t i, size, ti_is_be;
403 ti_is_be = bam_is_big_endian();
404 bgzf_write(fp, "TBI\1", 4);
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;
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
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);
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
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
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);
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);
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);
462 // write linear index (index2)
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
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);
478 static ti_index_t *ti_index_load_core(BGZF *fp)
483 ti_is_be = bam_is_big_endian();
485 fprintf(stderr, "[ti_index_load_core] fail to load index.\n");
488 bgzf_read(fp, magic, 4);
489 if (strncmp(magic, "TBI\1", 4)) {
490 fprintf(stderr, "[ti_index_load] wrong magic number.\n");
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));
500 bgzf_read(fp, &idx->conf, sizeof(ti_conf_t));
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);
509 { // read target names
514 bgzf_read(fp, &l, 4);
515 if (ti_is_be) bam_swap_endian_4p(&l);
517 bgzf_read(fp, buf, l);
518 str = calloc(1, sizeof(kstring_t));
519 for (i = j = 0; i < l; ++i) {
521 khint_t k = kh_put(s, idx->tname, strdup(str->s), &ret);
522 kh_value(idx->tname, k) = j++;
524 } else kputc(buf[i], str);
526 free(str->s); free(str); free(buf);
528 for (i = 0; i < idx->n; ++i) {
530 ti_lidx_t *index2 = idx->index2 + i;
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);
547 p->list = (pair64_t*)malloc(p->m * 16);
548 bgzf_read(fp, p->list, 16 * p->n);
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);
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);
564 for (j = 0; j < index2->n; ++j) bam_swap_endian_8p(&index2->offset[j]);
569 ti_index_t *ti_index_load_local(const char *fnidx)
572 fp = bgzf_open(fnidx, "r");
574 ti_index_t *idx = ti_index_load_core(fp);
581 static void download_from_remote(const char *url)
583 const int buf_size = 1 * 1024 * 1024;
589 if (strstr(url, "ftp://") != url && strstr(url, "http://") != url) return;
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");
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);
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);
609 knet_close(fp_remote);
612 static void download_from_remote(const char *url)
618 static char *get_local_version(const char *fn)
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)) {
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) {
633 fprintf(stderr, "[%s] downloading the index file...\n", __func__);
634 download_from_remote(url);
637 if (stat(fnidx, &sbuf) == 0) return fnidx;
638 free(fnidx); return 0;
641 const char **ti_seqname(const ti_index_t *idx, int *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);
653 ti_index_t *ti_index_load(const char *fn)
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);
664 int ti_index_build2(const char *fn, const ti_conf_t *conf, const char *_fnidx)
669 if ((fp = bgzf_open(fn, "r")) == 0) {
670 fprintf(stderr, "[ti_index_build2] fail to open the file: %s\n", fn);
673 idx = ti_index_core(fp, conf);
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");
681 fprintf(stderr, "[ti_index_build2] fail to create the index file.\n");
685 ti_index_save(idx, fpidx);
686 ti_index_destroy(idx);
692 int ti_index_build(const char *fn, const ti_conf_t *conf)
694 return ti_index_build2(fn, conf, 0);
697 /********************************************
698 * parse a region in the format chr:beg-end *
699 ********************************************/
701 int ti_get_tid(const ti_index_t *idx, const char *name)
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);
710 int ti_parse_region(const ti_index_t *idx, const char *str, int *tid, int *begin, int *end)
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];
720 for (i = 0; i != k; ++i) if (s[i] == ':') break;
722 if ((*tid = ti_get_tid(idx, s)) < 0) {
726 if (i == k) { /* dump the whole sequence */
727 *begin = 0; *end = 1<<29; free(s);
730 for (p = s + i + 1; i != k; ++i) if (s[i] == '-') break;
736 if (*begin > 0) --*begin;
738 if (*begin > *end) return -1;
742 /*******************************
743 * retrieve a specified region *
744 *******************************/
746 #define MAX_BIN 37450 // =(8^6-1)/7+1
748 static inline int reg2bins(uint32_t beg, uint32_t end, uint16_t list[MAX_BIN])
751 if (beg >= end) return 0;
752 if (end >= 1u<<29) end = 1u<<29;
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;
763 ti_iter_t ti_iter_first()
766 iter = calloc(1, sizeof(struct __ti_iter_t));
767 iter->from_first = 1;
771 ti_iter_t ti_iter_query(const ti_index_t *idx, int tid, int beg, int end)
774 int i, n_bins, n_off;
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;
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];
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;
806 free(bins); return iter;
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)) {
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];
818 free(bins); free(off); return iter;
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)
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];
840 iter->n_off = n_off; iter->off = off;
844 const char *ti_iter_read(BGZF *fp, ti_iter_t iter, int *len)
846 if (iter->finished) return 0;
847 if (iter->from_first) {
849 if ((ret = ti_readline(fp, &iter->str)) < 0) {
853 if (len) *len = iter->str.l;
857 if (iter->n_off == 0) return 0;
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);
869 if ((ret = ti_readline(fp, &iter->str)) >= 0) {
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;
879 } else break; // end of file
885 void ti_iter_destroy(ti_iter_t iter)
888 free(iter->str.s); free(iter->off);
893 int ti_fetch(BGZF *fp, const ti_index_t *idx, int tid, int beg, int end, void *data, ti_fetch_f func)
898 iter = ti_iter_query(idx, tid, beg, end);
899 while ((s = ti_iter_read(fp, iter, &len)) != 0)
901 ti_iter_destroy(iter);
909 tabix_t *ti_open(const char *fn, const char *fnidx)
913 if ((fp = bgzf_open(fn, "r")) == 0) return 0;
914 t = calloc(1, sizeof(tabix_t));
916 if (fnidx) t->fnidx = strdup(fnidx);
921 void ti_close(tabix_t *t)
925 if (t->idx) ti_index_destroy(t->idx);
926 free(t->fn); free(t->fnidx);
931 int ti_lazy_index_load(tabix_t *t)
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
941 ti_iter_t ti_queryi(tabix_t *t, int tid, int beg, int end)
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);
948 ti_iter_t ti_querys(tabix_t *t, const char *reg)
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);
957 ti_iter_t ti_query(tabix_t *t, const char *name, int beg, int end)
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);
967 const char *ti_read(tabix_t *t, ti_iter_t iter, int *len)
969 return ti_iter_read(t->fp, iter, len);