1 /// @file htslib/vcf.h
2 /// High-level VCF/BCF variant calling file operations.
3 /*
4     Copyright (C) 2012, 2013 Broad Institute.
5     Copyright (C) 2012-2014 Genome Research Ltd.
6 
7     Author: Heng Li <lh3@sanger.ac.uk>
8 
9 Permission is hereby granted, free of charge, to any person obtaining a copy
10 of this software and associated documentation files (the "Software"), to deal
11 in the Software without restriction, including without limitation the rights
12 to use, copy, modify, merge, publish, distribute, sublicense, and/or sell
13 copies of the Software, and to permit persons to whom the Software is
14 furnished to do so, subject to the following conditions:
15 
16 The above copyright notice and this permission notice shall be included in
17 all copies or substantial portions of the Software.
18 
19 THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
20 IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
21 FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
22 THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
23 LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
24 FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
25 DEALINGS IN THE SOFTWARE.  */
26 
27 /*
28     todo:
29         - make the function names consistent
30         - provide calls to abstract away structs as much as possible
31  */
32 
33 #ifndef HTSLIB_VCF_H
34 #define HTSLIB_VCF_H
35 
36 #include <stdint.h>
37 #include <limits.h>
38 #include <assert.h>
39 #include "hts.h"
40 #include "kstring.h"
41 #include "hts_defs.h"
42 #include "hts_endian.h"
43 
44 #ifdef __cplusplus
45 extern "C" {
46 #endif
47 
48 /*****************
49  * Header struct *
50  *****************/
51 
52 #define BCF_HL_FLT  0 // header line
53 #define BCF_HL_INFO 1
54 #define BCF_HL_FMT  2
55 #define BCF_HL_CTG  3
56 #define BCF_HL_STR  4 // structured header line TAG=<A=..,B=..>
57 #define BCF_HL_GEN  5 // generic header line
58 
59 #define BCF_HT_FLAG 0 // header type
60 #define BCF_HT_INT  1
61 #define BCF_HT_REAL 2
62 #define BCF_HT_STR  3
63 
64 #define BCF_VL_FIXED 0 // variable length
65 #define BCF_VL_VAR   1
66 #define BCF_VL_A     2
67 #define BCF_VL_G     3
68 #define BCF_VL_R     4
69 
70 /* === Dictionary ===
71 
72    The header keeps three dictonaries. The first keeps IDs in the
73    "FILTER/INFO/FORMAT" lines, the second keeps the sequence names and lengths
74    in the "contig" lines and the last keeps the sample names. bcf_hdr_t::dict[]
75    is the actual hash table, which is opaque to the end users. In the hash
76    table, the key is the ID or sample name as a C string and the value is a
77    bcf_idinfo_t struct. bcf_hdr_t::id[] points to key-value pairs in the hash
78    table in the order that they appear in the VCF header. bcf_hdr_t::n[] is the
79    size of the hash table or, equivalently, the length of the id[] arrays.
80 */
81 
82 #define BCF_DT_ID       0 // dictionary type
83 #define BCF_DT_CTG      1
84 #define BCF_DT_SAMPLE   2
85 
86 // Complete textual representation of a header line
87 typedef struct {
88     int type;       // One of the BCF_HL_* type
89     char *key;      // The part before '=', i.e. FILTER/INFO/FORMAT/contig/fileformat etc.
90     char *value;    // Set only for generic lines, NULL for FILTER/INFO, etc.
91     int nkeys;              // Number of structured fields
92     char **keys, **vals;    // The key=value pairs
93 } bcf_hrec_t;
94 
95 typedef struct {
96     uint32_t info[3];  // stores Number:20, var:4, Type:4, ColType:4 in info[0..2]
97                        // for BCF_HL_FLT,INFO,FMT and contig length in info[0] for BCF_HL_CTG
98     bcf_hrec_t *hrec[3];
99     int id;
100 } bcf_idinfo_t;
101 
102 typedef struct {
103     const char *key;
104     const bcf_idinfo_t *val;
105 } bcf_idpair_t;
106 
107 // Note that bcf_hdr_t structs must always be created via bcf_hdr_init()
108 typedef struct {
109     int32_t n[3];           // n:the size of the dictionary block in use, (allocated size, m, is below to preserve ABI)
110     bcf_idpair_t *id[3];
111     void *dict[3];          // ID dictionary, contig dict and sample dict
112     char **samples;
113     bcf_hrec_t **hrec;
114     int nhrec, dirty;
115     int ntransl, *transl[2];    // for bcf_translate()
116     int nsamples_ori;           // for bcf_hdr_set_samples()
117     uint8_t *keep_samples;
118     kstring_t mem;
119     int32_t m[3];          // m: allocated size of the dictionary block in use (see n above)
120 } bcf_hdr_t;
121 
122 extern uint8_t bcf_type_shift[];
123 
124 /**************
125  * VCF record *
126  **************/
127 
128 #define BCF_BT_NULL     0
129 #define BCF_BT_INT8     1
130 #define BCF_BT_INT16    2
131 #define BCF_BT_INT32    3
132 #define BCF_BT_FLOAT    5
133 #define BCF_BT_CHAR     7
134 
135 #define VCF_REF   0
136 #define VCF_SNP   1
137 #define VCF_MNP   2
138 #define VCF_INDEL 4
139 #define VCF_OTHER 8
140 #define VCF_BND   16    // breakend
141 
142 typedef struct {
143     int type, n;    // variant type and the number of bases affected, negative for deletions
144 } variant_t;
145 
146 typedef struct {
147     int id;             // id: numeric tag id, the corresponding string is bcf_hdr_t::id[BCF_DT_ID][$id].key
148     int n, size, type;  // n: number of values per-sample; size: number of bytes per-sample; type: one of BCF_BT_* types
149     uint8_t *p;         // same as vptr and vptr_* in bcf_info_t below
150     uint32_t p_len;
151     uint32_t p_off:31, p_free:1;
152 } bcf_fmt_t;
153 
154 typedef struct {
155     int key;        // key: numeric tag id, the corresponding string is bcf_hdr_t::id[BCF_DT_ID][$key].key
156     int type, len;  // type: one of BCF_BT_* types; len: vector length, 1 for scalars
157     union {
158         int32_t i; // integer value
159         float f;   // float value
160     } v1; // only set if $len==1; for easier access
161     uint8_t *vptr;          // pointer to data array in bcf1_t->shared.s, excluding the size+type and tag id bytes
162     uint32_t vptr_len;      // length of the vptr block or, when set, of the vptr_mod block, excluding offset
163     uint32_t vptr_off:31,   // vptr offset, i.e., the size of the INFO key plus size+type bytes
164             vptr_free:1;    // indicates that vptr-vptr_off must be freed; set only when modified and the new
165                             //    data block is bigger than the original
166 } bcf_info_t;
167 
168 
169 #define BCF1_DIRTY_ID  1
170 #define BCF1_DIRTY_ALS 2
171 #define BCF1_DIRTY_FLT 4
172 #define BCF1_DIRTY_INF 8
173 
174 typedef struct {
175     int m_fmt, m_info, m_id, m_als, m_allele, m_flt; // allocated size (high-water mark); do not change
176     int n_flt;  // Number of FILTER fields
177     int *flt;   // FILTER keys in the dictionary
178     char *id, *als;     // ID and REF+ALT block (\0-seperated)
179     char **allele;      // allele[0] is the REF (allele[] pointers to the als block); all null terminated
180     bcf_info_t *info;   // INFO
181     bcf_fmt_t *fmt;     // FORMAT and individual sample
182     variant_t *var;     // $var and $var_type set only when set_variant_types called
183     int n_var, var_type;
184     int shared_dirty;   // if set, shared.s must be recreated on BCF output
185     int indiv_dirty;    // if set, indiv.s must be recreated on BCF output
186 } bcf_dec_t;
187 
188 
189 #define BCF_ERR_CTG_UNDEF 1
190 #define BCF_ERR_TAG_UNDEF 2
191 #define BCF_ERR_NCOLS     4
192 #define BCF_ERR_LIMITS    8
193 #define BCF_ERR_CHAR     16
194 #define BCF_ERR_CTG_INVALID   32
195 #define BCF_ERR_TAG_INVALID   64
196 
197 /*
198     The bcf1_t structure corresponds to one VCF/BCF line. Reading from VCF file
199     is slower because the string is first to be parsed, packed into BCF line
200     (done in vcf_parse), then unpacked into internal bcf1_t structure. If it
201     is known in advance that some of the fields will not be required (notably
202     the sample columns), parsing of these can be skipped by setting max_unpack
203     appropriately.
204     Similarly, it is fast to output a BCF line because the columns (kept in
205     shared.s, indiv.s, etc.) are written directly by bcf_write, whereas a VCF
206     line must be formatted in vcf_format.
207  */
208 typedef struct {
209     int32_t rid;  // CHROM
210     int32_t pos;  // POS
211     int32_t rlen; // length of REF
212     float qual;   // QUAL
213     uint32_t n_info:16, n_allele:16;
214     uint32_t n_fmt:8, n_sample:24;
215     kstring_t shared, indiv;
216     bcf_dec_t d; // lazy evaluation: $d is not generated by bcf_read(), but by explicitly calling bcf_unpack()
217     int max_unpack;         // Set to BCF_UN_STR, BCF_UN_FLT, or BCF_UN_INFO to boost performance of vcf_parse when some of the fields won't be needed
218     int unpacked;           // remember what has been unpacked to allow calling bcf_unpack() repeatedly without redoing the work
219     int unpack_size[3];     // the original block size of ID, REF+ALT and FILTER
220     int errcode;    // one of BCF_ERR_* codes
221 } bcf1_t;
222 
223 /*******
224  * API *
225  *******/
226 
227     /***********************************************************************
228      *  BCF and VCF I/O
229      *
230      *  A note about naming conventions: htslib internally represents VCF
231      *  records as bcf1_t data structures, therefore most functions are
232      *  prefixed with bcf_. There are a few exceptions where the functions must
233      *  be aware of both BCF and VCF worlds, such as bcf_parse vs vcf_parse. In
234      *  these cases, functions prefixed with bcf_ are more general and work
235      *  with both BCF and VCF.
236      *
237      ***********************************************************************/
238 
239     /** These macros are defined only for consistency with other parts of htslib */
240     #define bcf_init1()         bcf_init()
241     #define bcf_read1(fp,h,v)   bcf_read((fp),(h),(v))
242     #define vcf_read1(fp,h,v)   vcf_read((fp),(h),(v))
243     #define bcf_write1(fp,h,v)  bcf_write((fp),(h),(v))
244     #define vcf_write1(fp,h,v)  vcf_write((fp),(h),(v))
245     #define bcf_destroy1(v)     bcf_destroy(v)
246     #define bcf_empty1(v)       bcf_empty(v)
247     #define vcf_parse1(s,h,v)   vcf_parse((s),(h),(v))
248     #define bcf_clear1(v)       bcf_clear(v)
249     #define vcf_format1(h,v,s)  vcf_format((h),(v),(s))
250 
251     /**
252      *  bcf_hdr_init() - create an empty BCF header.
253      *  @param mode    "r" or "w"
254      *
255      *  When opened for writing, the mandatory fileFormat and
256      *  FILTER=PASS lines are added automatically.
257      */
258     bcf_hdr_t *bcf_hdr_init(const char *mode);
259 
260     /** Destroy a BCF header struct */
261     void bcf_hdr_destroy(bcf_hdr_t *h);
262 
263     /** Initialize a bcf1_t object; equivalent to calloc(1, sizeof(bcf1_t)) */
264     bcf1_t *bcf_init(void);
265 
266     /** Deallocate a bcf1_t object */
267     void bcf_destroy(bcf1_t *v);
268 
269     /**
270      *  Same as bcf_destroy() but frees only the memory allocated by bcf1_t,
271      *  not the bcf1_t object itself.
272      */
273     void bcf_empty(bcf1_t *v);
274 
275     /**
276      *  Make the bcf1_t object ready for next read. Intended mostly for
277      *  internal use, the user should rarely need to call this function
278      *  directly.
279      */
280     void bcf_clear(bcf1_t *v);
281 
282 
283     /** bcf_open and vcf_open mode: please see hts_open() in hts.h */
284     typedef htsFile vcfFile;
285     #define bcf_open(fn, mode) hts_open((fn), (mode))
286     #define vcf_open(fn, mode) hts_open((fn), (mode))
287     #define bcf_close(fp) hts_close(fp)
288     #define vcf_close(fp) hts_close(fp)
289 
290     /** Reads VCF or BCF header */
291     bcf_hdr_t *bcf_hdr_read(htsFile *fp);
292 
293     /**
294      *  bcf_hdr_set_samples() - for more efficient VCF parsing when only one/few samples are needed
295      *  @samples: samples to include or exclude from file or as a comma-separated string.
296      *              LIST|FILE   .. select samples in list/file
297      *              ^LIST|FILE  .. exclude samples from list/file
298      *              -           .. include all samples
299      *              NULL        .. exclude all samples
300      *  @is_file: @samples is a file (1) or a comma-separated list (0)
301      *
302      *  The bottleneck of VCF reading is parsing of genotype fields. If the
303      *  reader knows in advance that only subset of samples is needed (possibly
304      *  no samples at all), the performance of bcf_read() can be significantly
305      *  improved by calling bcf_hdr_set_samples after bcf_hdr_read().
306      *  The function bcf_read() will subset the VCF/BCF records automatically
307      *  with the notable exception when reading records via bcf_itr_next().
308      *  In this case, bcf_subset_format() must be called explicitly, because
309      *  bcf_readrec() does not see the header.
310      *
311      *  Returns 0 on success, -1 on error or a positive integer if the list
312      *  contains samples not present in the VCF header. In such a case, the
313      *  return value is the index of the offending sample.
314      */
315     int bcf_hdr_set_samples(bcf_hdr_t *hdr, const char *samples, int is_file);
316     int bcf_subset_format(const bcf_hdr_t *hdr, bcf1_t *rec);
317 
318 
319     /** Writes VCF or BCF header */
320     int bcf_hdr_write(htsFile *fp, bcf_hdr_t *h);
321 
322     /**
323      * Parse VCF line contained in kstring and populate the bcf1_t struct
324      * The line must not end with \n or \r characters.
325      */
326     int vcf_parse(kstring_t *s, const bcf_hdr_t *h, bcf1_t *v);
327 
328     /** The opposite of vcf_parse. It should rarely be called directly, see vcf_write */
329     int vcf_format(const bcf_hdr_t *h, const bcf1_t *v, kstring_t *s);
330 
331     /**
332      *  bcf_read() - read next VCF or BCF record
333      *
334      *  Returns -1 on critical errors, 0 otherwise. On errors which are not
335      *  critical for reading, such as missing header definitions, v->errcode is
336      *  set to one of BCF_ERR* code and must be checked before calling
337      *  vcf_write().
338      */
339     int bcf_read(htsFile *fp, const bcf_hdr_t *h, bcf1_t *v);
340 
341     /**
342      *  bcf_unpack() - unpack/decode a BCF record (fills the bcf1_t::d field)
343      *
344      *  Note that bcf_unpack() must be called even when reading VCF. It is safe
345      *  to call the function repeatedly, it will not unpack the same field
346      *  twice.
347      */
348     #define BCF_UN_STR  1       // up to ALT inclusive
349     #define BCF_UN_FLT  2       // up to FILTER
350     #define BCF_UN_INFO 4       // up to INFO
351     #define BCF_UN_SHR  (BCF_UN_STR|BCF_UN_FLT|BCF_UN_INFO) // all shared information
352     #define BCF_UN_FMT  8                           // unpack format and each sample
353     #define BCF_UN_IND  BCF_UN_FMT                  // a synonymo of BCF_UN_FMT
354     #define BCF_UN_ALL  (BCF_UN_SHR|BCF_UN_FMT)     // everything
355     int bcf_unpack(bcf1_t *b, int which);
356 
357     /*
358      *  bcf_dup() - create a copy of BCF record.
359      *
360      *  Note that bcf_unpack() must be called on the returned copy as if it was
361      *  obtained from bcf_read(). Also note that bcf_dup() calls bcf_sync1(src)
362      *  internally to reflect any changes made by bcf_update_* functions.
363      */
364     bcf1_t *bcf_dup(bcf1_t *src);
365     bcf1_t *bcf_copy(bcf1_t *dst, bcf1_t *src);
366 
367     /**
368      *  bcf_write() - write one VCF or BCF record. The type is determined at the open() call.
369      */
370     int bcf_write(htsFile *fp, bcf_hdr_t *h, bcf1_t *v);
371 
372     /**
373      *  The following functions work only with VCFs and should rarely be called
374      *  directly. Usually one wants to use their bcf_* alternatives, which work
375      *  transparently with both VCFs and BCFs.
376      */
377     bcf_hdr_t *vcf_hdr_read(htsFile *fp);
378     int vcf_hdr_write(htsFile *fp, const bcf_hdr_t *h);
379     int vcf_read(htsFile *fp, const bcf_hdr_t *h, bcf1_t *v);
380     int vcf_write(htsFile *fp, const bcf_hdr_t *h, bcf1_t *v);
381 
382     /** Helper function for the bcf_itr_next() macro; internal use, ignore it */
383     int bcf_readrec(BGZF *fp, void *null, void *v, int *tid, int *beg, int *end);
384 
385 
386 
387     /**************************************************************************
388      *  Header querying and manipulation routines
389      **************************************************************************/
390 
391     /** Create a new header using the supplied template */
392     bcf_hdr_t *bcf_hdr_dup(const bcf_hdr_t *hdr);
393 
394     /**
395      *  Copy header lines from src to dst if not already present in dst. See also bcf_translate().
396      *  Returns 0 on success or sets a bit on error:
397      *      1 .. conflicting definitions of tag length
398      *      // todo
399      */
400     int bcf_hdr_combine(bcf_hdr_t *dst, const bcf_hdr_t *src) HTS_DEPRECATED("Please use bcf_hdr_merge instead");
401 
402     /**
403      *  bcf_hdr_merge() - copy header lines from src to dst, see also bcf_translate()
404      *  @param dst: the destination header to be merged into, NULL on the first pass
405      *  @param src: the source header
406      *
407      *  Notes:
408      *      - use as:
409      *          bcf_hdr_t *dst = NULL;
410      *          for (i=0; i<nsrc; i++) dst = bcf_hdr_merge(dst,src[i]);
411      *
412      *      - bcf_hdr_merge() replaces bcf_hdr_combine() which had a problem when
413      *      combining multiple BCF headers. The current bcf_hdr_combine()
414      *      does not have this problem, but became slow when used for many files.
415      */
416     bcf_hdr_t *bcf_hdr_merge(bcf_hdr_t *dst, const bcf_hdr_t *src);
417 
418     /**
419      *  bcf_hdr_add_sample() - add a new sample.
420      *  @param sample:  sample name to be added
421      */
422     int bcf_hdr_add_sample(bcf_hdr_t *hdr, const char *sample);
423 
424     /** Read VCF header from a file and update the header */
425     int bcf_hdr_set(bcf_hdr_t *hdr, const char *fname);
426 
427     /// Appends formatted header text to _str_.
428     /** If _is_bcf_ is zero, `IDX` fields are discarded.
429      *  @return 0 if successful, or negative if an error occurred
430      *  @since 1.4
431      */
432     int bcf_hdr_format(const bcf_hdr_t *hdr, int is_bcf, kstring_t *str);
433 
434     /** Returns formatted header (newly allocated string) and its length,
435      *  excluding the terminating \0. If is_bcf parameter is unset, IDX
436      *  fields are discarded.
437      *  @deprecated Use bcf_hdr_format() instead as it can handle huge headers.
438      */
439     char *bcf_hdr_fmt_text(const bcf_hdr_t *hdr, int is_bcf, int *len)
440         HTS_DEPRECATED("use bcf_hdr_format() instead");
441 
442     /** Append new VCF header line, returns 0 on success */
443     int bcf_hdr_append(bcf_hdr_t *h, const char *line);
444     int bcf_hdr_printf(bcf_hdr_t *h, const char *format, ...);
445 
446     /** VCF version, e.g. VCFv4.2 */
447     const char *bcf_hdr_get_version(const bcf_hdr_t *hdr);
448     void bcf_hdr_set_version(bcf_hdr_t *hdr, const char *version);
449 
450     /**
451      *  bcf_hdr_remove() - remove VCF header tag
452      *  @param type:      one of BCF_HL_*
453      *  @param key:       tag name or NULL to remove all tags of the given type
454      */
455     void bcf_hdr_remove(bcf_hdr_t *h, int type, const char *key);
456 
457     /**
458      *  bcf_hdr_subset() - creates a new copy of the header removing unwanted samples
459      *  @param n:        number of samples to keep
460      *  @param samples:  names of the samples to keep
461      *  @param imap:     mapping from index in @samples to the sample index in the original file
462      *
463      *  Sample names not present in h0 are ignored. The number of unmatched samples can be checked
464      *  by comparing n and bcf_hdr_nsamples(out_hdr).
465      *  This function can be used to reorder samples.
466      *  See also bcf_subset() which subsets individual records.
467      */
468     bcf_hdr_t *bcf_hdr_subset(const bcf_hdr_t *h0, int n, char *const* samples, int *imap);
469 
470     /** Creates a list of sequence names. It is up to the caller to free the list (but not the sequence names) */
471     const char **bcf_hdr_seqnames(const bcf_hdr_t *h, int *nseqs);
472 
473     /** Get number of samples */
474     #define bcf_hdr_nsamples(hdr) (hdr)->n[BCF_DT_SAMPLE]
475 
476 
477     /** The following functions are for internal use and should rarely be called directly */
478     int bcf_hdr_parse(bcf_hdr_t *hdr, char *htxt);
479     int bcf_hdr_sync(bcf_hdr_t *h);
480     bcf_hrec_t *bcf_hdr_parse_line(const bcf_hdr_t *h, const char *line, int *len);
481     void bcf_hrec_format(const bcf_hrec_t *hrec, kstring_t *str);
482     int bcf_hdr_add_hrec(bcf_hdr_t *hdr, bcf_hrec_t *hrec);
483     /**
484      *  bcf_hdr_get_hrec() - get header line info
485      *  @param type:  one of the BCF_HL_* types: FLT,INFO,FMT,CTG,STR,GEN
486      *  @param key:   the header key for generic lines (e.g. "fileformat"), any field
487      *                  for structured lines, typically "ID".
488      *  @param value: the value which pairs with key. Can be be NULL for BCF_HL_GEN
489      *  @param str_class: the class of BCF_HL_STR line (e.g. "ALT" or "SAMPLE"), otherwise NULL
490      */
491     bcf_hrec_t *bcf_hdr_get_hrec(const bcf_hdr_t *hdr, int type, const char *key, const char *value, const char *str_class);
492     bcf_hrec_t *bcf_hrec_dup(bcf_hrec_t *hrec);
493     void bcf_hrec_add_key(bcf_hrec_t *hrec, const char *str, int len);
494     void bcf_hrec_set_val(bcf_hrec_t *hrec, int i, const char *str, int len, int is_quoted);
495     int bcf_hrec_find_key(bcf_hrec_t *hrec, const char *key);
496     void hrec_add_idx(bcf_hrec_t *hrec, int idx);
497     void bcf_hrec_destroy(bcf_hrec_t *hrec);
498 
499 
500 
501     /**************************************************************************
502      *  Individual record querying and manipulation routines
503      **************************************************************************/
504 
505     /** See the description of bcf_hdr_subset() */
506     int bcf_subset(const bcf_hdr_t *h, bcf1_t *v, int n, int *imap);
507 
508     /**
509      *  bcf_translate() - translate tags ids to be consistent with different header. This function
510      *                    is useful when lines from multiple VCF need to be combined.
511      *  @dst_hdr:   the destination header, to be used in bcf_write(), see also bcf_hdr_combine()
512      *  @src_hdr:   the source header, used in bcf_read()
513      *  @src_line:  line obtained by bcf_read()
514      */
515     int bcf_translate(const bcf_hdr_t *dst_hdr, bcf_hdr_t *src_hdr, bcf1_t *src_line);
516 
517     /**
518      *  bcf_get_variant_type[s]()  - returns one of VCF_REF, VCF_SNP, etc
519      */
520     int bcf_get_variant_types(bcf1_t *rec);
521     int bcf_get_variant_type(bcf1_t *rec, int ith_allele);
522     int bcf_is_snp(bcf1_t *v);
523 
524     /**
525      *  bcf_update_filter() - sets the FILTER column
526      *  @flt_ids:  The filter IDs to set, numeric IDs returned by bcf_hdr_id2int(hdr, BCF_DT_ID, "PASS")
527      *  @n:        Number of filters. If n==0, all filters are removed
528      */
529     int bcf_update_filter(const bcf_hdr_t *hdr, bcf1_t *line, int *flt_ids, int n);
530     /**
531      *  bcf_add_filter() - adds to the FILTER column
532      *  @flt_id:   filter ID to add, numeric ID returned by bcf_hdr_id2int(hdr, BCF_DT_ID, "PASS")
533      *
534      *  If flt_id is PASS, all existing filters are removed first. If other than PASS, existing PASS is removed.
535      */
536     int bcf_add_filter(const bcf_hdr_t *hdr, bcf1_t *line, int flt_id);
537     /**
538      *  bcf_remove_filter() - removes from the FILTER column
539      *  @flt_id:   filter ID to remove, numeric ID returned by bcf_hdr_id2int(hdr, BCF_DT_ID, "PASS")
540      *  @pass:     when set to 1 and no filters are present, set to PASS
541      */
542     int bcf_remove_filter(const bcf_hdr_t *hdr, bcf1_t *line, int flt_id, int pass);
543     /**
544      *  Returns 1 if present, 0 if absent, or -1 if filter does not exist. "PASS" and "." can be used interchangeably.
545      */
546     int bcf_has_filter(const bcf_hdr_t *hdr, bcf1_t *line, char *filter);
547     /**
548      *  bcf_update_alleles() and bcf_update_alleles_str() - update REF and ALLT column
549      *  @alleles:           Array of alleles
550      *  @nals:              Number of alleles
551      *  @alleles_string:    Comma-separated alleles, starting with the REF allele
552      */
553     int bcf_update_alleles(const bcf_hdr_t *hdr, bcf1_t *line, const char **alleles, int nals);
554     int bcf_update_alleles_str(const bcf_hdr_t *hdr, bcf1_t *line, const char *alleles_string);
555 
556     /**
557       *  bcf_update_id() - sets new ID string
558       *  bcf_add_id() - adds to the ID string checking for duplicates
559       */
560     int bcf_update_id(const bcf_hdr_t *hdr, bcf1_t *line, const char *id);
561     int bcf_add_id(const bcf_hdr_t *hdr, bcf1_t *line, const char *id);
562 
563     /*
564      *  bcf_update_info_*() - functions for updating INFO fields
565      *  @hdr:       the BCF header
566      *  @line:      VCF line to be edited
567      *  @key:       the INFO tag to be updated
568      *  @values:    pointer to the array of values. Pass NULL to remove the tag.
569      *  @n:         number of values in the array. When set to 0, the INFO tag is removed
570      *
571      *  The @string in bcf_update_info_flag() is optional, @n indicates whether
572      *  the flag is set or removed.
573      *
574      *  Returns 0 on success or negative value on error.
575      */
576     #define bcf_update_info_int32(hdr,line,key,values,n)   bcf_update_info((hdr),(line),(key),(values),(n),BCF_HT_INT)
577     #define bcf_update_info_float(hdr,line,key,values,n)   bcf_update_info((hdr),(line),(key),(values),(n),BCF_HT_REAL)
578     #define bcf_update_info_flag(hdr,line,key,string,n)    bcf_update_info((hdr),(line),(key),(string),(n),BCF_HT_FLAG)
579     #define bcf_update_info_string(hdr,line,key,string)    bcf_update_info((hdr),(line),(key),(string),1,BCF_HT_STR)
580     int bcf_update_info(const bcf_hdr_t *hdr, bcf1_t *line, const char *key, const void *values, int n, int type);
581 
582     /*
583      *  bcf_update_format_*() - functions for updating FORMAT fields
584      *  @values:    pointer to the array of values, the same number of elements
585      *              is expected for each sample. Missing values must be padded
586      *              with bcf_*_missing or bcf_*_vector_end values.
587      *  @n:         number of values in the array. If n==0, existing tag is removed.
588      *
589      *  The function bcf_update_format_string() is a higher-level (slower) variant of
590      *  bcf_update_format_char(). The former accepts array of \0-terminated strings
591      *  whereas the latter requires that the strings are collapsed into a single array
592      *  of fixed-length strings. In case of strings with variable length, shorter strings
593      *  can be \0-padded. Note that the collapsed strings passed to bcf_update_format_char()
594      *  are not \0-terminated.
595      *
596      *  Returns 0 on success or negative value on error.
597      */
598     #define bcf_update_format_int32(hdr,line,key,values,n) bcf_update_format((hdr),(line),(key),(values),(n),BCF_HT_INT)
599     #define bcf_update_format_float(hdr,line,key,values,n) bcf_update_format((hdr),(line),(key),(values),(n),BCF_HT_REAL)
600     #define bcf_update_format_char(hdr,line,key,values,n) bcf_update_format((hdr),(line),(key),(values),(n),BCF_HT_STR)
601     #define bcf_update_genotypes(hdr,line,gts,n) bcf_update_format((hdr),(line),"GT",(gts),(n),BCF_HT_INT)     // See bcf_gt_ macros below
602     int bcf_update_format_string(const bcf_hdr_t *hdr, bcf1_t *line, const char *key, const char **values, int n);
603     int bcf_update_format(const bcf_hdr_t *hdr, bcf1_t *line, const char *key, const void *values, int n, int type);
604 
605     // Macros for setting genotypes correctly, for use with bcf_update_genotypes only; idx corresponds
606     // to VCF's GT (1-based index to ALT or 0 for the reference allele) and val is the opposite, obtained
607     // from bcf_get_genotypes() below.
608     #define bcf_gt_phased(idx)      (((idx)+1)<<1|1)
609     #define bcf_gt_unphased(idx)    (((idx)+1)<<1)
610     #define bcf_gt_missing          0
611     #define bcf_gt_is_missing(val)  ((val)>>1 ? 0 : 1)
612     #define bcf_gt_is_phased(idx)   ((idx)&1)
613     #define bcf_gt_allele(val)      (((val)>>1)-1)
614 
615     /** Conversion between alleles indexes to Number=G genotype index (assuming diploid, all 0-based) */
616     #define bcf_alleles2gt(a,b) ((a)>(b)?((a)*((a)+1)/2+(b)):((b)*((b)+1)/2+(a)))
bcf_gt2alleles(int igt,int * a,int * b)617     static inline void bcf_gt2alleles(int igt, int *a, int *b)
618     {
619         int k = 0, dk = 1;
620         while ( k<igt ) { dk++; k += dk; }
621         *b = dk - 1; *a = igt - k + *b;
622     }
623 
624     /**
625      * bcf_get_fmt() - returns pointer to FORMAT's field data
626      * @header: for access to BCF_DT_ID dictionary
627      * @line:   VCF line obtained from vcf_parse1
628      * @fmt:    one of GT,PL,...
629      *
630      * Returns bcf_fmt_t* if the call succeeded, or returns NULL when the field
631      * is not available.
632      */
633     bcf_fmt_t *bcf_get_fmt(const bcf_hdr_t *hdr, bcf1_t *line, const char *key);
634     bcf_info_t *bcf_get_info(const bcf_hdr_t *hdr, bcf1_t *line, const char *key);
635 
636     /**
637      * bcf_get_*_id() - returns pointer to FORMAT/INFO field data given the header index instead of the string ID
638      * @line: VCF line obtained from vcf_parse1
639      * @id:  The header index for the tag, obtained from bcf_hdr_id2int()
640      *
641      * Returns bcf_fmt_t* / bcf_info_t*. These functions do not check if the index is valid
642      * as their goal is to avoid the header lookup.
643      */
644     bcf_fmt_t *bcf_get_fmt_id(bcf1_t *line, const int id);
645     bcf_info_t *bcf_get_info_id(bcf1_t *line, const int id);
646 
647     /**
648      *  bcf_get_info_*() - get INFO values, integers or floats
649      *  @hdr:       BCF header
650      *  @line:      BCF record
651      *  @tag:       INFO tag to retrieve
652      *  @dst:       *dst is pointer to a memory location, can point to NULL
653      *  @ndst:      pointer to the size of allocated memory
654      *
655      *  Returns negative value on error or the number of written values
656      *  (including missing values) on success. bcf_get_info_string() returns
657      *  on success the number of characters written excluding the null-
658      *  terminating byte. bcf_get_info_flag() returns 1 when flag is set or 0
659      *  if not.
660      *
661      *  List of return codes:
662      *      -1 .. no such INFO tag defined in the header
663      *      -2 .. clash between types defined in the header and encountered in the VCF record
664      *      -3 .. tag is not present in the VCF record
665      */
666     #define bcf_get_info_int32(hdr,line,tag,dst,ndst)  bcf_get_info_values(hdr,line,tag,(void**)(dst),ndst,BCF_HT_INT)
667     #define bcf_get_info_float(hdr,line,tag,dst,ndst)  bcf_get_info_values(hdr,line,tag,(void**)(dst),ndst,BCF_HT_REAL)
668     #define bcf_get_info_string(hdr,line,tag,dst,ndst) bcf_get_info_values(hdr,line,tag,(void**)(dst),ndst,BCF_HT_STR)
669     #define bcf_get_info_flag(hdr,line,tag,dst,ndst)   bcf_get_info_values(hdr,line,tag,(void**)(dst),ndst,BCF_HT_FLAG)
670     int bcf_get_info_values(const bcf_hdr_t *hdr, bcf1_t *line, const char *tag, void **dst, int *ndst, int type);
671 
672     /**
673      *  bcf_get_format_*() - same as bcf_get_info*() above
674      *
675      *  The function bcf_get_format_string() is a higher-level (slower) variant of bcf_get_format_char().
676      *  see the description of bcf_update_format_string() and bcf_update_format_char() above.
677      *  Unlike other bcf_get_format__*() functions, bcf_get_format_string() allocates two arrays:
678      *  a single block of \0-terminated strings collapsed into a single array and an array of pointers
679      *  to these strings. Both arrays must be cleaned by the user.
680      *
681      *  Returns negative value on error or the number of written values on success.
682      *
683      *  Use the returned number of written values for accessing valid entries of dst, as ndst is only a
684      *  watermark that can be higher than the returned value, i.e. the end of dst can contain carry-over
685      *  values from previous calls to bcf_get_format_*() on lines with more values per sample.
686      *
687      *  Example:
688      *      int ndst = 0; char **dst = NULL;
689      *      if ( bcf_get_format_string(hdr, line, "XX", &dst, &ndst) > 0 )
690      *          for (i=0; i<bcf_hdr_nsamples(hdr); i++) printf("%s\n", dst[i]);
691      *      free(dst[0]); free(dst);
692      *
693      *  Example:
694      *      int i, j, ngt, nsmpl = bcf_hdr_nsamples(hdr);
695      *      int32_t *gt_arr = NULL, ngt_arr = 0;
696      *
697      *      ngt = bcf_get_genotypes(hdr, line, &gt_arr, &ngt_arr);
698      *      if ( ngt<=0 ) return; // GT not present
699      *
700      *      int max_ploidy = ngt/nsmpl;
701      *      for (i=0; i<nsmpl; i++)
702      *      {
703      *        int32_t *ptr = gt + i*max_ploidy;
704      *        for (j=0; j<max_ploidy; j++)
705      *        {
706      *           // if true, the sample has smaller ploidy
707      *           if ( ptr[j]==bcf_int32_vector_end ) break;
708      *
709      *           // missing allele
710      *           if ( bcf_gt_is_missing(ptr[j]) ) continue;
711      *
712      *           // the VCF 0-based allele index
713      *           int allele_index = bcf_gt_allele(ptr[j]);
714      *
715      *           // is phased?
716      *           int is_phased = bcf_gt_is_phased(ptr[j]);
717      *
718      *           // .. do something ..
719      *         }
720      *      }
721      *      free(gt_arr);
722      *
723      */
724     #define bcf_get_format_int32(hdr,line,tag,dst,ndst)  bcf_get_format_values(hdr,line,tag,(void**)(dst),ndst,BCF_HT_INT)
725     #define bcf_get_format_float(hdr,line,tag,dst,ndst)  bcf_get_format_values(hdr,line,tag,(void**)(dst),ndst,BCF_HT_REAL)
726     #define bcf_get_format_char(hdr,line,tag,dst,ndst)   bcf_get_format_values(hdr,line,tag,(void**)(dst),ndst,BCF_HT_STR)
727     #define bcf_get_genotypes(hdr,line,dst,ndst)         bcf_get_format_values(hdr,line,"GT",(void**)(dst),ndst,BCF_HT_INT)
728     int bcf_get_format_string(const bcf_hdr_t *hdr, bcf1_t *line, const char *tag, char ***dst, int *ndst);
729     int bcf_get_format_values(const bcf_hdr_t *hdr, bcf1_t *line, const char *tag, void **dst, int *ndst, int type);
730 
731 
732 
733     /**************************************************************************
734      *  Helper functions
735      **************************************************************************/
736 
737     /**
738      *  bcf_hdr_id2int() - Translates string into numeric ID
739      *  bcf_hdr_int2id() - Translates numeric ID into string
740      *  @type:     one of BCF_DT_ID, BCF_DT_CTG, BCF_DT_SAMPLE
741      *  @id:       tag name, such as: PL, DP, GT, etc.
742      *
743      *  Returns -1 if string is not in dictionary, otherwise numeric ID which identifies
744      *  fields in BCF records.
745      */
746     int bcf_hdr_id2int(const bcf_hdr_t *hdr, int type, const char *id);
747     #define bcf_hdr_int2id(hdr,type,int_id) ((hdr)->id[type][int_id].key)
748 
749     /**
750      *  bcf_hdr_name2id() - Translates sequence names (chromosomes) into numeric ID
751      *  bcf_hdr_id2name() - Translates numeric ID to sequence name
752      */
bcf_hdr_name2id(const bcf_hdr_t * hdr,const char * id)753     static inline int bcf_hdr_name2id(const bcf_hdr_t *hdr, const char *id) { return bcf_hdr_id2int(hdr, BCF_DT_CTG, id); }
bcf_hdr_id2name(const bcf_hdr_t * hdr,int rid)754     static inline const char *bcf_hdr_id2name(const bcf_hdr_t *hdr, int rid) { return hdr->id[BCF_DT_CTG][rid].key; }
bcf_seqname(const bcf_hdr_t * hdr,bcf1_t * rec)755     static inline const char *bcf_seqname(const bcf_hdr_t *hdr, bcf1_t *rec) { return hdr->id[BCF_DT_CTG][rec->rid].key; }
756 
757     /**
758      *  bcf_hdr_id2*() - Macros for accessing bcf_idinfo_t
759      *  @type:      one of BCF_HL_FLT, BCF_HL_INFO, BCF_HL_FMT
760      *  @int_id:    return value of bcf_hdr_id2int, must be >=0
761      *
762      *  The returned values are:
763      *     bcf_hdr_id2length   ..  whether the number of values is fixed or variable, one of BCF_VL_*
764      *     bcf_hdr_id2number   ..  the number of values, 0xfffff for variable length fields
765      *     bcf_hdr_id2type     ..  the field type, one of BCF_HT_*
766      *     bcf_hdr_id2coltype  ..  the column type, one of BCF_HL_*
767      *
768      *  Notes: Prior to using the macros, the presence of the info should be
769      *  tested with bcf_hdr_idinfo_exists().
770      */
771     #define bcf_hdr_id2length(hdr,type,int_id)  ((hdr)->id[BCF_DT_ID][int_id].val->info[type]>>8 & 0xf)
772     #define bcf_hdr_id2number(hdr,type,int_id)  ((hdr)->id[BCF_DT_ID][int_id].val->info[type]>>12)
773     #define bcf_hdr_id2type(hdr,type,int_id)    ((hdr)->id[BCF_DT_ID][int_id].val->info[type]>>4 & 0xf)
774     #define bcf_hdr_id2coltype(hdr,type,int_id) ((hdr)->id[BCF_DT_ID][int_id].val->info[type] & 0xf)
775     #define bcf_hdr_idinfo_exists(hdr,type,int_id)  ((int_id<0 || bcf_hdr_id2coltype(hdr,type,int_id)==0xf) ? 0 : 1)
776     #define bcf_hdr_id2hrec(hdr,dict_type,col_type,int_id)    ((hdr)->id[(dict_type)==BCF_DT_CTG?BCF_DT_CTG:BCF_DT_ID][int_id].val->hrec[(dict_type)==BCF_DT_CTG?0:(col_type)])
777 
778     void bcf_fmt_array(kstring_t *s, int n, int type, void *data);
779     uint8_t *bcf_fmt_sized_array(kstring_t *s, uint8_t *ptr);
780 
781     void bcf_enc_vchar(kstring_t *s, int l, const char *a);
782     void bcf_enc_vint(kstring_t *s, int n, int32_t *a, int wsize);
783     void bcf_enc_vfloat(kstring_t *s, int n, float *a);
784 
785 
786     /**************************************************************************
787      *  BCF index
788      *
789      *  Note that these functions work with BCFs only. See synced_bcf_reader.h
790      *  which provides (amongst other things) an API to work transparently with
791      *  both indexed BCFs and VCFs.
792      **************************************************************************/
793 
794     #define bcf_itr_destroy(iter) hts_itr_destroy(iter)
795     #define bcf_itr_queryi(idx, tid, beg, end) hts_itr_query((idx), (tid), (beg), (end), bcf_readrec)
796     #define bcf_itr_querys(idx, hdr, s) hts_itr_querys((idx), (s), (hts_name2id_f)(bcf_hdr_name2id), (hdr), hts_itr_query, bcf_readrec)
797     #define bcf_itr_next(htsfp, itr, r) hts_itr_next((htsfp)->fp.bgzf, (itr), (r), 0)
798     #define bcf_index_load(fn) hts_idx_load(fn, HTS_FMT_CSI)
799     #define bcf_index_seqnames(idx, hdr, nptr) hts_idx_seqnames((idx),(nptr),(hts_id2name_f)(bcf_hdr_id2name),(hdr))
800 
801     hts_idx_t *bcf_index_load2(const char *fn, const char *fnidx);
802 
803     /**
804      *  bcf_index_build() - Generate and save an index file
805      *  @fn:         Input VCF/BCF filename
806      *  @min_shift:  Positive to generate CSI, or 0 to generate TBI
807      *
808      *  Returns 0 if successful, or negative if an error occurred.
809      *
810      *  List of error codes:
811      *      -1 .. indexing failed
812      *      -2 .. opening @fn failed
813      *      -3 .. format not indexable
814      *      -4 .. failed to create and/or save the index
815      */
816     int bcf_index_build(const char *fn, int min_shift);
817 
818     /**
819      *  bcf_index_build2() - Generate and save an index to a specific file
820      *  @fn:         Input VCF/BCF filename
821      *  @fnidx:      Output filename, or NULL to add .csi/.tbi to @fn
822      *  @min_shift:  Positive to generate CSI, or 0 to generate TBI
823      *
824      *  Returns 0 if successful, or negative if an error occurred.
825      *
826      *  List of error codes:
827      *      -1 .. indexing failed
828      *      -2 .. opening @fn failed
829      *      -3 .. format not indexable
830      *      -4 .. failed to create and/or save the index
831      */
832     int bcf_index_build2(const char *fn, const char *fnidx, int min_shift);
833 
834     /**
835      *  bcf_index_build3() - Generate and save an index to a specific file
836      *  @fn:         Input VCF/BCF filename
837      *  @fnidx:      Output filename, or NULL to add .csi/.tbi to @fn
838      *  @min_shift:  Positive to generate CSI, or 0 to generate TBI
839      *  @n_threads:  Number of VCF/BCF decoder threads
840      *
841      *  Returns 0 if successful, or negative if an error occurred.
842      *
843      *  List of error codes:
844      *      -1 .. indexing failed
845      *      -2 .. opening @fn failed
846      *      -3 .. format not indexable
847      *      -4 .. failed to create and/or save the index
848      */
849      int bcf_index_build3(const char *fn, const char *fnidx, int min_shift, int n_threads);
850 
851 /*******************
852  * Typed value I/O *
853  *******************/
854 
855 /*
856     Note that in contrast with BCFv2.1 specification, HTSlib implementation
857     allows missing values in vectors. For integer types, the values 0x80,
858     0x8000, 0x80000000 are interpreted as missing values and 0x81, 0x8001,
859     0x80000001 as end-of-vector indicators.  Similarly for floats, the value of
860     0x7F800001 is interpreted as a missing value and 0x7F800002 as an
861     end-of-vector indicator.
862     Note that the end-of-vector byte is not part of the vector.
863 
864     This trial BCF version (v2.2) is compatible with the VCF specification and
865     enables to handle correctly vectors with different ploidy in presence of
866     missing values.
867  */
868 #define bcf_int8_vector_end  (INT8_MIN+1)
869 #define bcf_int16_vector_end (INT16_MIN+1)
870 #define bcf_int32_vector_end (INT32_MIN+1)
871 #define bcf_str_vector_end   0
872 #define bcf_int8_missing     INT8_MIN
873 #define bcf_int16_missing    INT16_MIN
874 #define bcf_int32_missing    INT32_MIN
875 #define bcf_str_missing      0x07
876 extern uint32_t bcf_float_vector_end;
877 extern uint32_t bcf_float_missing;
bcf_float_set(float * ptr,uint32_t value)878 static inline void bcf_float_set(float *ptr, uint32_t value)
879 {
880     union { uint32_t i; float f; } u;
881     u.i = value;
882     *ptr = u.f;
883 }
884 #define bcf_float_set_vector_end(x) bcf_float_set(&(x),bcf_float_vector_end)
885 #define bcf_float_set_missing(x)    bcf_float_set(&(x),bcf_float_missing)
bcf_float_is_missing(float f)886 static inline int bcf_float_is_missing(float f)
887 {
888     union { uint32_t i; float f; } u;
889     u.f = f;
890     return u.i==bcf_float_missing ? 1 : 0;
891 }
bcf_float_is_vector_end(float f)892 static inline int bcf_float_is_vector_end(float f)
893 {
894     union { uint32_t i; float f; } u;
895     u.f = f;
896     return u.i==bcf_float_vector_end ? 1 : 0;
897 }
898 
bcf_format_gt(bcf_fmt_t * fmt,int isample,kstring_t * str)899 static inline void bcf_format_gt(bcf_fmt_t *fmt, int isample, kstring_t *str)
900 {
901     #define BRANCH(type_t, missing, vector_end) { \
902         type_t *ptr = (type_t*) (fmt->p + isample*fmt->size); \
903         int i; \
904         for (i=0; i<fmt->n && ptr[i]!=vector_end; i++) \
905         { \
906             if ( i ) kputc("/|"[ptr[i]&1], str); \
907             if ( !(ptr[i]>>1) ) kputc('.', str); \
908             else kputw((ptr[i]>>1) - 1, str); \
909         } \
910         if (i == 0) kputc('.', str); \
911     }
912     switch (fmt->type) {
913         case BCF_BT_INT8:  BRANCH(int8_t,  bcf_int8_missing, bcf_int8_vector_end); break;
914         case BCF_BT_INT16: BRANCH(int16_t, bcf_int16_missing, bcf_int16_vector_end); break;
915         case BCF_BT_INT32: BRANCH(int32_t, bcf_int32_missing, bcf_int32_vector_end); break;
916         case BCF_BT_NULL:  kputc('.', str); break;
917         default: fprintf(stderr,"FIXME: type %d in bcf_format_gt?\n", fmt->type); abort(); break;
918     }
919     #undef BRANCH
920 }
921 
bcf_enc_size(kstring_t * s,int size,int type)922 static inline void bcf_enc_size(kstring_t *s, int size, int type)
923 {
924     if (size >= 15) {
925         kputc(15<<4|type, s);
926         if (size >= 128) {
927             if (size >= 32768) {
928                 int32_t x = size;
929                 kputc(1<<4|BCF_BT_INT32, s);
930                 kputsn((char*)&x, 4, s);
931             } else {
932                 int16_t x = size;
933                 kputc(1<<4|BCF_BT_INT16, s);
934                 kputsn((char*)&x, 2, s);
935             }
936         } else {
937             kputc(1<<4|BCF_BT_INT8, s);
938             kputc(size, s);
939         }
940     } else kputc(size<<4|type, s);
941 }
942 
bcf_enc_inttype(long x)943 static inline int bcf_enc_inttype(long x)
944 {
945     if (x <= INT8_MAX && x > bcf_int8_missing) return BCF_BT_INT8;
946     if (x <= INT16_MAX && x > bcf_int16_missing) return BCF_BT_INT16;
947     return BCF_BT_INT32;
948 }
949 
bcf_enc_int1(kstring_t * s,int32_t x)950 static inline void bcf_enc_int1(kstring_t *s, int32_t x)
951 {
952     if (x == bcf_int32_vector_end) {
953         bcf_enc_size(s, 1, BCF_BT_INT8);
954         kputc(bcf_int8_vector_end, s);
955     } else if (x == bcf_int32_missing) {
956         bcf_enc_size(s, 1, BCF_BT_INT8);
957         kputc(bcf_int8_missing, s);
958     } else if (x <= INT8_MAX && x > bcf_int8_missing) {
959         bcf_enc_size(s, 1, BCF_BT_INT8);
960         kputc(x, s);
961     } else if (x <= INT16_MAX && x > bcf_int16_missing) {
962         int16_t z = x;
963         bcf_enc_size(s, 1, BCF_BT_INT16);
964         kputsn((char*)&z, 2, s);
965     } else {
966         int32_t z = x;
967         bcf_enc_size(s, 1, BCF_BT_INT32);
968         kputsn((char*)&z, 4, s);
969     }
970 }
971 
bcf_dec_int1(const uint8_t * p,int type,uint8_t ** q)972 static inline int32_t bcf_dec_int1(const uint8_t *p, int type, uint8_t **q)
973 {
974     if (type == BCF_BT_INT8) {
975         *q = (uint8_t*)p + 1;
976         return le_to_i8(p);
977     } else if (type == BCF_BT_INT16) {
978         *q = (uint8_t*)p + 2;
979         return le_to_i16(p);
980     } else {
981         *q = (uint8_t*)p + 4;
982         return le_to_i32(p);
983     }
984 }
985 
bcf_dec_typed_int1(const uint8_t * p,uint8_t ** q)986 static inline int32_t bcf_dec_typed_int1(const uint8_t *p, uint8_t **q)
987 {
988     return bcf_dec_int1(p + 1, *p&0xf, q);
989 }
990 
bcf_dec_size(const uint8_t * p,uint8_t ** q,int * type)991 static inline int32_t bcf_dec_size(const uint8_t *p, uint8_t **q, int *type)
992 {
993     *type = *p & 0xf;
994     if (*p>>4 != 15) {
995         *q = (uint8_t*)p + 1;
996         return *p>>4;
997     } else return bcf_dec_typed_int1(p + 1, q);
998 }
999 
1000 #ifdef __cplusplus
1001 }
1002 #endif
1003 
1004 #endif
1005