@@ -35,6 +35,7 @@ DEALINGS IN THE SOFTWARE. */
3535#include "htslib/faidx.h"
3636#include "htslib/hfile.h"
3737#include "htslib/khash.h"
38+ #include "htslib/kstring.h"
3839
3940typedef struct {
4041 int32_t line_len , line_blen ;
@@ -80,9 +81,8 @@ static inline void fai_insert_index(faidx_t *idx, const char *name, int len, int
8081
8182faidx_t * fai_build_core (BGZF * bgzf )
8283{
83- char * name ;
84+ kstring_t name = { 0 , 0 , NULL } ;
8485 int c ;
85- int l_name , m_name ;
8686 int line_len , line_blen , state ;
8787 int l1 , l2 ;
8888 faidx_t * idx ;
@@ -91,7 +91,6 @@ faidx_t *fai_build_core(BGZF *bgzf)
9191
9292 idx = (faidx_t * )calloc (1 , sizeof (faidx_t ));
9393 idx -> hash = kh_init (s );
94- name = 0 ; l_name = m_name = 0 ;
9594 len = line_len = line_blen = -1 ; state = 0 ; l1 = l2 = -1 ; offset = 0 ;
9695 while ( (c = bgzf_getc (bgzf ))>=0 ) {
9796 if (c == '\n' ) { // an empty line
@@ -103,30 +102,25 @@ faidx_t *fai_build_core(BGZF *bgzf)
103102 }
104103 if (c == '>' ) { // fasta header
105104 if (len >= 0 )
106- fai_insert_index (idx , name , len , line_len , line_blen , offset );
107- l_name = 0 ;
108- while ( (c = bgzf_getc (bgzf ))>=0 && !isspace (c )) {
109- if (m_name < l_name + 2 ) {
110- m_name = l_name + 2 ;
111- kroundup32 (m_name );
112- name = (char * )realloc (name , m_name );
113- }
114- name [l_name ++ ] = c ;
115- }
116- name [l_name ] = '\0' ;
105+ fai_insert_index (idx , name .s , len , line_len , line_blen , offset );
106+
107+ name .l = 0 ;
108+ while ((c = bgzf_getc (bgzf )) >= 0 )
109+ if (! isspace (c )) kputc_ (c , & name );
110+ else if (name .l > 0 || c == '\n' ) break ;
111+ kputsn ("" , 0 , & name );
112+
117113 if ( c < 0 ) {
118114 fprintf (stderr , "[fai_build_core] the last entry has no sequence\n" );
119- free (name ); fai_destroy (idx );
120- return 0 ;
115+ goto fail ;
121116 }
122117 if (c != '\n' ) while ( (c = bgzf_getc (bgzf ))>=0 && c != '\n' );
123118 state = 1 ; len = 0 ;
124119 offset = bgzf_utell (bgzf );
125120 } else {
126121 if (state == 3 ) {
127- fprintf (stderr , "[fai_build_core] inlined empty line is not allowed in sequence '%s'.\n" , name );
128- free (name ); fai_destroy (idx );
129- return 0 ;
122+ fprintf (stderr , "[fai_build_core] inlined empty line is not allowed in sequence '%s'.\n" , name .s );
123+ goto fail ;
130124 }
131125 if (state == 2 ) state = 3 ;
132126 l1 = l2 = 0 ;
@@ -135,9 +129,8 @@ faidx_t *fai_build_core(BGZF *bgzf)
135129 if (isgraph (c )) ++ l2 ;
136130 } while ( (c = bgzf_getc (bgzf ))>=0 && c != '\n' );
137131 if (state == 3 && l2 ) {
138- fprintf (stderr , "[fai_build_core] different line length in sequence '%s'.\n" , name );
139- free (name ); fai_destroy (idx );
140- return 0 ;
132+ fprintf (stderr , "[fai_build_core] different line length in sequence '%s'.\n" , name .s );
133+ goto fail ;
141134 }
142135 ++ l1 ; len += l2 ;
143136 if (state == 1 ) line_len = l1 , line_blen = l2 , state = 0 ;
@@ -146,15 +139,19 @@ faidx_t *fai_build_core(BGZF *bgzf)
146139 }
147140 }
148141 }
149- if ( name )
150- fai_insert_index (idx , name , len , line_len , line_blen , offset );
142+
143+ if (len >= 0 )
144+ fai_insert_index (idx , name .s , len , line_len , line_blen , offset );
151145 else
152- {
153- free (idx );
154- return NULL ;
155- }
156- free (name );
146+ goto fail ;
147+
148+ free (name .s );
157149 return idx ;
150+
151+ fail :
152+ free (name .s );
153+ fai_destroy (idx );
154+ return NULL ;
158155}
159156
160157void fai_save (const faidx_t * fai , FILE * fp )
@@ -229,6 +226,7 @@ int fai_build(const char *fn)
229226 if ( !fai )
230227 {
231228 if ( bgzf -> is_compressed && bgzf -> is_gzip ) fprintf (stderr ,"Cannot index files compressed with gzip, please use bgzip\n" );
229+ bgzf_close (bgzf );
232230 free (str );
233231 return -1 ;
234232 }
0 commit comments