1 /// @file htslib/synced_bcf_reader.h 2 /// Stream through multiple VCF files. 3 /* 4 Copyright (C) 2012-2014 Genome Research Ltd. 5 6 Author: Petr Danecek <pd3@sanger.ac.uk> 7 8 Permission is hereby granted, free of charge, to any person obtaining a copy 9 of this software and associated documentation files (the "Software"), to deal 10 in the Software without restriction, including without limitation the rights 11 to use, copy, modify, merge, publish, distribute, sublicense, and/or sell 12 copies of the Software, and to permit persons to whom the Software is 13 furnished to do so, subject to the following conditions: 14 15 The above copyright notice and this permission notice shall be included in 16 all copies or substantial portions of the Software. 17 18 THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR 19 IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY, 20 FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL 21 THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER 22 LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING 23 FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER 24 DEALINGS IN THE SOFTWARE. */ 25 26 /* 27 The synced_bcf_reader allows to keep multiple VCFs open and stream them 28 using the next_line iterator in a seamless matter without worrying about 29 chromosomes and synchronizing the sites. This is used by vcfcheck to 30 compare multiple VCFs simultaneously and is used also for merging, 31 creating intersections, etc. 32 33 The synced_bcf_reader also provides API for reading indexed BCF/VCF, 34 hiding differences in BCF/VCF opening, indexing and reading. 35 36 37 Example of usage: 38 39 bcf_srs_t *sr = bcf_sr_init(); 40 bcf_sr_set_opt(sr, BCF_SR_PAIR_LOGIC, BCF_SR_PAIR_BOTH_REF); 41 bcf_sr_set_opt(sr, BCF_SR_REQUIRE_IDX); 42 for (i=0; i<nfiles; i++) 43 bcf_sr_add_reader(sr,files[i]); 44 while ( bcf_sr_next_line(sr) ) 45 { 46 for (i=0; i<nfiles; i++) 47 { 48 bcf1_t *line = bcf_sr_get_line(sr,i); 49 ... 50 } 51 } 52 if ( sr->errnum ) error("Error: %s\n", bcf_sr_strerror(sr->errnum)); 53 bcf_sr_destroy(sr); 54 */ 55 56 #ifndef HTSLIB_SYNCED_BCF_READER_H 57 #define HTSLIB_SYNCED_BCF_READER_H 58 59 #include "hts.h" 60 #include "vcf.h" 61 #include "tbx.h" 62 63 #ifdef __cplusplus 64 extern "C" { 65 #endif 66 67 /* 68 When reading multiple files in paralel, duplicate records within each 69 file will be reordered and offered in intuitive order. For example, 70 when reading two files, each with unsorted SNP and indel record, the 71 reader should return the SNP records together and the indel records 72 together. The logic of compatible records can vary depending on the 73 application and can be set using the PAIR_* defined below. 74 75 The COLLAPSE_* definitions will be deprecated in future versions, please 76 use the PAIR_* definitions instead. 77 */ 78 #define COLLAPSE_NONE 0 // require the exact same set of alleles in all files 79 #define COLLAPSE_SNPS 1 // allow different alleles, as long as they all are SNPs 80 #define COLLAPSE_INDELS 2 // the same as above, but with indels 81 #define COLLAPSE_ANY 4 // any combination of alleles can be returned by bcf_sr_next_line() 82 #define COLLAPSE_SOME 8 // at least some of the ALTs must match 83 #define COLLAPSE_BOTH (COLLAPSE_SNPS|COLLAPSE_INDELS) 84 85 #define BCF_SR_PAIR_SNPS (1<<0) // allow different alleles, as long as they all are SNPs 86 #define BCF_SR_PAIR_INDELS (1<<1) // the same as above, but with indels 87 #define BCF_SR_PAIR_ANY (1<<2) // any combination of alleles can be returned by bcf_sr_next_line() 88 #define BCF_SR_PAIR_SOME (1<<3) // at least some of multiallelic ALTs must match. Implied by all the others with the exception of EXACT 89 #define BCF_SR_PAIR_SNP_REF (1<<4) // allow REF-only records with SNPs 90 #define BCF_SR_PAIR_INDEL_REF (1<<5) // allow REF-only records with indels 91 #define BCF_SR_PAIR_EXACT (1<<6) // require the exact same set of alleles in all files 92 #define BCF_SR_PAIR_BOTH (BCF_SR_PAIR_SNPS|BCF_SR_PAIR_INDELS) 93 #define BCF_SR_PAIR_BOTH_REF (BCF_SR_PAIR_SNPS|BCF_SR_PAIR_INDELS|BCF_SR_PAIR_SNP_REF|BCF_SR_PAIR_INDEL_REF) 94 95 typedef enum 96 { 97 BCF_SR_REQUIRE_IDX, 98 BCF_SR_PAIR_LOGIC // combination of the PAIR_* values above 99 } 100 bcf_sr_opt_t; 101 102 typedef struct _bcf_sr_regions_t 103 { 104 // for reading from tabix-indexed file (big data) 105 tbx_t *tbx; // tabix index 106 hts_itr_t *itr; // tabix iterator 107 kstring_t line; // holder of the current line, set only when reading from tabix-indexed files 108 htsFile *file; 109 char *fname; 110 int is_bin; // is open in binary mode (tabix access) 111 char **als; // parsed alleles if targets_als set and _regions_match_alleles called 112 kstring_t als_str; // block of parsed alleles 113 int nals, mals; // number of set alleles and the size of allocated array 114 int als_type; // alleles type, currently VCF_SNP or VCF_INDEL 115 116 // user handler to deal with skipped regions without a counterpart in VCFs 117 void (*missed_reg_handler)(struct _bcf_sr_regions_t *, void *); 118 void *missed_reg_data; 119 120 // for in-memory regions (small data) 121 struct _region_t *regs; // the regions 122 123 // shared by both tabix-index and in-memory regions 124 void *seq_hash; // keys: sequence names, values: index to seqs 125 char **seq_names; // sequence names 126 int nseqs; // number of sequences (chromosomes) in the file 127 int iseq; // current position: chr name, index to snames 128 int start, end; // current position: start, end of the region (0-based) 129 int prev_seq, prev_start; 130 } 131 bcf_sr_regions_t; 132 133 typedef struct 134 { 135 htsFile *file; 136 tbx_t *tbx_idx; 137 hts_idx_t *bcf_idx; 138 bcf_hdr_t *header; 139 hts_itr_t *itr; 140 char *fname; 141 bcf1_t **buffer; // cached VCF records. First is the current record synced across the reader 142 int nbuffer, mbuffer; // number of cached records (including the current record); number of allocated records 143 int nfilter_ids, *filter_ids; // -1 for ".", otherwise filter id as returned by bcf_hdr_id2int 144 int *samples, n_smpl; // list of columns in the order consistent with bcf_srs_t.samples 145 } 146 bcf_sr_t; 147 148 typedef enum 149 { 150 open_failed, not_bgzf, idx_load_failed, file_type_error, api_usage_error, 151 header_error, no_eof, no_memory, vcf_parse_error, bcf_read_error 152 } 153 bcf_sr_error; 154 155 typedef struct 156 { 157 // Parameters controlling the logic 158 int collapse; // Do not access directly, use bcf_sr_set_pairing_logic() instead 159 char *apply_filters; // If set, sites where none of the FILTER strings is listed 160 // will be skipped. Active only at the time of 161 // initialization, that is during the add_reader() 162 // calls. Therefore, each reader can be initialized with different 163 // filters. 164 int require_index; // Some tools do not need random access 165 int max_unpack; // When reading VCFs and knowing some fields will not be needed, boost performance of vcf_parse1 166 int *has_line; // Corresponds to return value of bcf_sr_next_line but is not limited by sizeof(int). Use bcf_sr_has_line macro to query. 167 bcf_sr_error errnum; 168 169 // Auxiliary data 170 bcf_sr_t *readers; 171 int nreaders; 172 int streaming; // reading mode: index-jumping or streaming 173 int explicit_regs; // was the list of regions se by bcf_sr_set_regions or guessed from tabix index? 174 char **samples; // List of samples 175 bcf_sr_regions_t *regions, *targets; // see bcf_sr_set_[targets|regions] for description 176 int targets_als; // subset to targets not only by position but also by alleles? 177 int targets_exclude; 178 kstring_t tmps; 179 int n_smpl; 180 181 int n_threads; // Simple multi-threaded decoding / encoding. 182 htsThreadPool *p; // Our pool, but it can be used by others if needed. 183 void *aux; // Opaque auxiliary data 184 } 185 bcf_srs_t; 186 187 /** Init bcf_srs_t struct */ 188 bcf_srs_t *bcf_sr_init(void); 189 190 /** Destroy bcf_srs_t struct */ 191 void bcf_sr_destroy(bcf_srs_t *readers); 192 193 char *bcf_sr_strerror(int errnum); 194 195 int bcf_sr_set_opt(bcf_srs_t *readers, bcf_sr_opt_t opt, ...); 196 197 198 /** 199 * bcf_sr_set_threads() - allocates a thread-pool for use by the synced reader. 200 * @n_threads: size of thread pool 201 * 202 * Returns 0 if the call succeeded, or <0 on error. 203 */ 204 int bcf_sr_set_threads(bcf_srs_t *files, int n_threads); 205 206 /** Deallocates thread memory, if owned by us. */ 207 void bcf_sr_destroy_threads(bcf_srs_t *files); 208 209 /** 210 * bcf_sr_add_reader() - open new reader 211 * @readers: holder of the open readers 212 * @fname: the VCF file 213 * 214 * Returns 1 if the call succeeded, or 0 on error. 215 * 216 * See also the bcf_srs_t data structure for parameters controlling 217 * the reader's logic. 218 */ 219 int bcf_sr_add_reader(bcf_srs_t *readers, const char *fname); 220 void bcf_sr_remove_reader(bcf_srs_t *files, int i); 221 222 /** 223 * bcf_sr_next_line() - the iterator 224 * @readers: holder of the open readers 225 * 226 * Returns the number of readers which have the current line 227 * (bcf_sr_t.buffer[0]) set at this position. Use the bcf_sr_has_line macro to 228 * determine which of the readers are set. 229 */ 230 int bcf_sr_next_line(bcf_srs_t *readers); 231 #define bcf_sr_has_line(readers, i) (readers)->has_line[i] 232 #define bcf_sr_get_line(_readers, i) ((_readers)->has_line[i] ? ((_readers)->readers[i].buffer[0]) : NULL) 233 #define bcf_sr_swap_line(_readers, i, lieu) { bcf1_t *tmp = lieu; lieu = (_readers)->readers[i].buffer[0]; (_readers)->readers[i].buffer[0] = tmp; } 234 #define bcf_sr_region_done(_readers,i) (!(_readers)->has_line[i] && !(_readers)->readers[i].nbuffer ? 1 : 0) 235 #define bcf_sr_get_header(_readers, i) (_readers)->readers[i].header 236 #define bcf_sr_get_reader(_readers, i) &((_readers)->readers[i]) 237 238 239 /** 240 * bcf_sr_seek() - set all readers to selected position 241 * @seq: sequence name; NULL to seek to start 242 * @pos: 0-based coordinate 243 */ 244 int bcf_sr_seek(bcf_srs_t *readers, const char *seq, int pos); 245 246 /** 247 * bcf_sr_set_samples() - sets active samples 248 * @readers: holder of the open readers 249 * @samples: this can be one of: file name with one sample per line; 250 * or column-separated list of samples; or '-' for a list of 251 * samples shared by all files. If first character is the 252 * exclamation mark, all but the listed samples are included. 253 * @is_file: 0: list of samples; 1: file with sample names 254 * 255 * Returns 1 if the call succeeded, or 0 on error. 256 */ 257 int bcf_sr_set_samples(bcf_srs_t *readers, const char *samples, int is_file); 258 259 /** 260 * bcf_sr_set_targets(), bcf_sr_set_regions() - init targets/regions 261 * @readers: holder of the open readers 262 * @targets: list of regions, one-based and inclusive. 263 * @is_fname: 0: targets is a comma-separated list of regions (chr,chr:from-to) 264 * 1: targets is a tabix indexed file with a list of regions 265 * (<chr,pos> or <chr,from,to>) 266 * 267 * Returns 0 if the call succeeded, or -1 on error. 268 * 269 * Both functions behave the same way, unlisted positions will be skipped by 270 * bcf_sr_next_line(). However, there is an important difference: regions use 271 * index to jump to desired positions while targets streams the whole files 272 * and merely skip unlisted positions. 273 * 274 * Moreover, bcf_sr_set_targets() accepts an optional parameter $alleles which 275 * is intepreted as a 1-based column index in the tab-delimited file where 276 * alleles are listed. This in principle enables to perform the COLLAPSE_* 277 * logic also with tab-delimited files. However, the current implementation 278 * considers the alleles merely as a suggestion for prioritizing one of possibly 279 * duplicate VCF lines. It is up to the caller to examine targets->als if 280 * perfect match is sought after. Note that the duplicate positions in targets 281 * file are currently not supported. 282 * Targets (but not regions) can be prefixed with "^" to request logical complement, 283 * for example "^X,Y,MT" indicates that sequences X, Y and MT should be skipped. 284 */ 285 int bcf_sr_set_targets(bcf_srs_t *readers, const char *targets, int is_file, int alleles); 286 int bcf_sr_set_regions(bcf_srs_t *readers, const char *regions, int is_file); 287 288 289 290 /* 291 * bcf_sr_regions_init() 292 * @regions: regions can be either a comma-separated list of regions 293 * (chr|chr:pos|chr:from-to|chr:from-) or VCF, BED, or 294 * tab-delimited file (the default). Uncompressed files 295 * are stored in memory while bgzip-compressed and tabix-indexed 296 * region files are streamed. 297 * @is_file: 0: regions is a comma-separated list of regions 298 * (chr|chr:pos|chr:from-to|chr:from-) 299 * 1: VCF, BED or tab-delimited file 300 * @chr, from, to: 301 * Column indexes of chromosome, start position and end position 302 * in the tab-delimited file. The positions are 1-based and 303 * inclusive. 304 * These parameters are ignored when reading from VCF, BED or 305 * tabix-indexed files. When end position column is not present, 306 * supply 'from' in place of 'to'. When 'to' is negative, first 307 * abs(to) will be attempted and if that fails, 'from' will be used 308 * instead. 309 */ 310 bcf_sr_regions_t *bcf_sr_regions_init(const char *regions, int is_file, int chr, int from, int to); 311 void bcf_sr_regions_destroy(bcf_sr_regions_t *regions); 312 313 /* 314 * bcf_sr_regions_seek() - seek to the chromosome block 315 * 316 * Returns 0 on success or -1 on failure. Sets reg->seq appropriately and 317 * reg->start,reg->end to -1. 318 */ 319 int bcf_sr_regions_seek(bcf_sr_regions_t *regions, const char *chr); 320 321 /* 322 * bcf_sr_regions_next() - retrieves next region. Returns 0 on success and -1 323 * when all regions have been read. The fields reg->seq, reg->start and 324 * reg->end are filled with the genomic coordinates on succes or with 325 * NULL,-1,-1 when no region is available. The coordinates are 0-based, 326 * inclusive. 327 */ 328 int bcf_sr_regions_next(bcf_sr_regions_t *reg); 329 330 /* 331 * bcf_sr_regions_overlap() - checks if the interval <start,end> overlaps any of 332 * the regions, the coordinates are 0-based, inclusive. The coordinate queries 333 * must come in ascending order. 334 * 335 * Returns 0 if the position is in regions; -1 if the position is not in the 336 * regions and more regions exist; -2 if not in the regions and there are no more 337 * regions left. 338 */ 339 int bcf_sr_regions_overlap(bcf_sr_regions_t *reg, const char *seq, int start, int end); 340 341 /* 342 * bcf_sr_regions_flush() - calls repeatedly regs->missed_reg_handler() until 343 * all remaining records are processed. 344 */ 345 void bcf_sr_regions_flush(bcf_sr_regions_t *regs); 346 347 #ifdef __cplusplus 348 } 349 #endif 350 351 #endif 352