static void _prune_sites(vcfbuf_t *buf, int flush_all) { int nbuf = flush_all ? buf->rbuf.n : buf->rbuf.n - 1; if ( nbuf > buf->prune.mvrec ) { buf->prune.idx = (int*) realloc(buf->prune.idx, nbuf*sizeof(int)); buf->prune.vrec = (vcfrec_t**) realloc(buf->prune.vrec, nbuf*sizeof(vcfrec_t*)); buf->prune.mvrec = nbuf; } // set allele frequency and prepare buffer for sorting int i,k,irec = 0; for (i=-1; rbuf_next(&buf->rbuf,&i) && irec<nbuf; ) { bcf1_t *line = buf->vcf[i].rec; if ( line->n_allele > buf->prune.mac ) { buf->prune.ac = (int*) realloc(buf->prune.ac, line->n_allele*sizeof(*buf->prune.ac)); buf->prune.mac = line->n_allele; } if ( !buf->vcf[i].af_set ) { buf->vcf[i].af = 0; if ( buf->prune.af_tag ) { if ( bcf_get_info_float(buf->hdr,line,buf->prune.af_tag,&buf->prune.farr, &buf->prune.mfarr) > 0 ) buf->vcf[i].af = buf->prune.farr[0]; } else if ( bcf_calc_ac(buf->hdr, line, buf->prune.ac, BCF_UN_INFO|BCF_UN_FMT) ) { int ntot = buf->prune.ac[0], nalt = 0; for (k=1; k<line->n_allele; k++) nalt += buf->prune.ac[k]; buf->vcf[i].af = ntot ? (float)nalt/ntot : 0; } buf->vcf[i].af_set = 1; } buf->vcf[i].idx = irec; buf->prune.vrec[irec++] = &buf->vcf[i]; } // sort by allele frequency, low AF will be removed preferentially qsort(buf->prune.vrec, nbuf, sizeof(*buf->prune.vrec), cmpvrec); // sort the rbuf indexes to be pruned descendently so that j-th rbuf index // is removed before i-th index if i<j int nprune = nbuf - buf->prune.max_sites; for (i=0; i<nprune; i++) buf->prune.idx[i] = buf->prune.vrec[i]->idx; qsort(buf->prune.idx, nprune, sizeof(int), cmpint_desc); for (i=0; i<nprune; i++) rbuf_remove_kth(&buf->rbuf, vcfrec_t, buf->prune.idx[i], buf->vcf); }
int process(bcf1_t *rec) { hts_expand(int,rec->n_allele,marr,arr); int ret = bcf_calc_ac(in_hdr,rec,arr,BCF_UN_FMT); if ( ret>0 ) { int i, an = 0; for (i=0; i<rec->n_allele; i++) an += arr[i]; bcf_update_info_int32(out_hdr, rec, "AN", &an, 1); bcf_update_info_int32(out_hdr, rec, "AC", arr+1, rec->n_allele-1); } return 0; }
static void buffered_filters(args_t *args, bcf1_t *line) { /** * The logic of SnpGap=3. The SNPs at positions 1 and 7 are filtered, * positions 0 and 8 are not: * 0123456789 * ref .G.GT..G.. * del .A.G-..A.. * Here the positions 1 and 6 are filtered, 0 and 7 are not: * 0123-456789 * ref .G.G-..G.. * ins .A.GT..A.. * * The logic of IndelGap=2. The second indel is filtered: * 012345678901 * ref .GT.GT..GT.. * del .G-.G-..G-.. * And similarly here, the second is filtered: * 01 23 456 78 * ref .A-.A-..A-.. * ins .AT.AT..AT.. */ // To avoid additional data structure, we abuse bcf1_t's var and var_type records. const int SnpGap_set = VCF_OTHER<<1; const int IndelGap_set = VCF_OTHER<<2; const int IndelGap_flush = VCF_OTHER<<3; int var_type = 0, i; if ( line ) { // Still on the same chromosome? int ilast = rbuf_last(&args->rbuf); if ( ilast>=0 && line->rid != args->rbuf_lines[ilast]->rid ) flush_buffer(args, args->rbuf.n); // new chromosome, flush everything rbuf_expand0(&args->rbuf,bcf1_t*,args->rbuf.n,args->rbuf_lines); // Insert the new record in the buffer. The line would be overwritten in // the next bcf_sr_next_line call, therefore we need to swap it with an // unused one ilast = rbuf_append(&args->rbuf); if ( !args->rbuf_lines[ilast] ) args->rbuf_lines[ilast] = bcf_init1(); SWAP(bcf1_t*, args->files->readers[0].buffer[0], args->rbuf_lines[ilast]); var_type = bcf_get_variant_types(line); // Find out the size of an indel. The indel boundaries are based on REF // (POS+1,POS+rlen-1). This is not entirely correct: mpileup likes to // output REF=CAGAGAGAGA, ALT=CAGAGAGAGAGA where REF=C,ALT=CGA could be // used. This filter is therefore more strict and may remove some valid // SNPs. int len = 1; if ( var_type & VCF_INDEL ) { for (i=1; i<line->n_allele; i++) if ( len < 1-line->d.var[i].n ) len = 1-line->d.var[i].n; } // Set the REF allele's length to max deletion length or to 1 if a SNP or an insertion. line->d.var[0].n = len; } int k_flush = 1; if ( args->indel_gap ) { k_flush = 0; // Find indels which are too close to each other int last_to = -1; for (i=-1; rbuf_next(&args->rbuf,&i); ) { bcf1_t *rec = args->rbuf_lines[i]; int rec_from = rec->pos; if ( last_to!=-1 && last_to < rec_from ) break; k_flush++; if ( !(rec->d.var_type & VCF_INDEL) ) continue; rec->d.var_type |= IndelGap_set; last_to = args->indel_gap + rec->pos + rec->d.var[0].n - 1; } if ( i==args->rbuf.f && line && last_to!=-1 ) k_flush = 0; if ( k_flush || !line ) { // Select the best indel from the cluster of k_flush indels int k = 0, max_ac = -1, imax_ac = -1; for (i=-1; rbuf_next(&args->rbuf,&i) && k<k_flush; ) { k++; bcf1_t *rec = args->rbuf_lines[i]; if ( !(rec->d.var_type & IndelGap_set) ) continue; hts_expand(int, rec->n_allele, args->ntmpi, args->tmpi); int ret = bcf_calc_ac(args->hdr, rec, args->tmpi, BCF_UN_ALL); if ( imax_ac==-1 || (ret && max_ac < args->tmpi[1]) ) { max_ac = args->tmpi[1]; imax_ac = i; } } // Filter all but the best indel (with max AF or first if AF not available) k = 0; for (i=-1; rbuf_next(&args->rbuf,&i) && k<k_flush; ) { k++; bcf1_t *rec = args->rbuf_lines[i]; if ( !(rec->d.var_type & IndelGap_set) ) continue; rec->d.var_type |= IndelGap_flush; if ( i!=imax_ac ) bcf_add_filter(args->hdr, args->rbuf_lines[i], args->IndelGap_id); } }
int subset_vcf(args_t *args, bcf1_t *line) { if ( args->min_alleles && line->n_allele < args->min_alleles ) return 0; // min alleles if ( args->max_alleles && line->n_allele > args->max_alleles ) return 0; // max alleles if (args->novel || args->known) { if ( args->novel && (line->d.id[0]!='.' || line->d.id[1]!=0) ) return 0; // skip sites which are known, ID != '.' if ( args->known && line->d.id[0]=='.' && line->d.id[1]==0 ) return 0; // skip sites which are novel, ID == '.' } if (args->include || args->exclude) { int line_type = bcf_get_variant_types(line); if ( args->include && !(line_type&args->include) ) return 0; // include only given variant types if ( args->exclude && line_type&args->exclude ) return 0; // exclude given variant types } if ( args->filter ) { int ret = filter_test(args->filter, line, NULL); if ( args->filter_logic==FLT_INCLUDE ) { if ( !ret ) return 0; } else if ( ret ) return 0; } hts_expand(int, line->n_allele, args->mac, args->ac); int i, an = 0, non_ref_ac = 0; if (args->calc_ac) { bcf_calc_ac(args->hdr, line, args->ac, BCF_UN_INFO|BCF_UN_FMT); // get original AC and AN values from INFO field if available, otherwise calculate for (i=1; i<line->n_allele; i++) non_ref_ac += args->ac[i]; for (i=0; i<line->n_allele; i++) an += args->ac[i]; } if (args->n_samples) { int non_ref_ac_sub = 0, *ac_sub = (int*) calloc(line->n_allele,sizeof(int)); bcf_subset(args->hdr, line, args->n_samples, args->imap); if (args->calc_ac) { bcf_calc_ac(args->hsub, line, ac_sub, BCF_UN_FMT); // recalculate AC and AN an = 0; for (i=0; i<line->n_allele; i++) { args->ac[i] = ac_sub[i]; an += ac_sub[i]; } for (i=1; i<line->n_allele; i++) non_ref_ac_sub += ac_sub[i]; if (args->private_vars) { if (args->private_vars == FLT_INCLUDE && !(non_ref_ac_sub > 0 && non_ref_ac == non_ref_ac_sub)) { free(ac_sub); return 0; } // select private sites if (args->private_vars == FLT_EXCLUDE && non_ref_ac_sub > 0 && non_ref_ac == non_ref_ac_sub) { free(ac_sub); return 0; } // exclude private sites } non_ref_ac = non_ref_ac_sub; } free(ac_sub); } bcf_fmt_t *gt_fmt; if ( args->gt_type && (gt_fmt=bcf_get_fmt(args->hdr,line,"GT")) ) { int nhet = 0, nhom = 0, nmiss = 0; for (i=0; i<bcf_hdr_nsamples(args->hdr); i++) { int type = bcf_gt_type(gt_fmt,i,NULL,NULL); if ( type==GT_HET_RA || type==GT_HET_AA ) { if ( args->gt_type==GT_NO_HET ) return 0; nhet = 1; } else if ( type==GT_UNKN ) { if ( args->gt_type==GT_NO_MISSING ) return 0; nmiss = 1; } else { if ( args->gt_type==GT_NO_HOM ) return 0; nhom = 1; } } if ( args->gt_type==GT_NEED_HOM && !nhom ) return 0; else if ( args->gt_type==GT_NEED_HET && !nhet ) return 0; else if ( args->gt_type==GT_NEED_MISSING && !nmiss ) return 0; } int minor_ac = 0; int major_ac = 0; if ( args->calc_ac ) { minor_ac = args->ac[0]; major_ac = args->ac[0]; for (i=1; i<line->n_allele; i++){ if (args->ac[i] < minor_ac) { minor_ac = args->ac[i]; } if (args->ac[i] > major_ac) { major_ac = args->ac[i]; } } } if (args->min_ac) { if (args->min_ac_type == ALLELE_NONREF && args->min_ac>non_ref_ac) return 0; // min AC else if (args->min_ac_type == ALLELE_MINOR && args->min_ac>minor_ac) return 0; // min minor AC else if (args->min_ac_type == ALLELE_ALT1 && args->min_ac>args->ac[1]) return 0; // min 1st alternate AC else if (args->min_ac_type == ALLELE_MAJOR && args->min_ac > major_ac) return 0; // min major AC else if (args->min_ac_type == ALLELE_NONMAJOR && args->min_ac > an-major_ac) return 0; // min non-major AC } if (args->max_ac) { if (args->max_ac_type == ALLELE_NONREF && args->max_ac<non_ref_ac) return 0; // max AC else if (args->max_ac_type == ALLELE_MINOR && args->max_ac<minor_ac) return 0; // max minor AC else if (args->max_ac_type == ALLELE_ALT1 && args->max_ac<args->ac[1]) return 0; // max 1st alternate AC else if (args->max_ac_type == ALLELE_MAJOR && args->max_ac < major_ac) return 0; // max major AC else if (args->max_ac_type == ALLELE_NONMAJOR && args->max_ac < an-major_ac) return 0; // max non-major AC } if (args->min_af) { if (an == 0) return 0; // freq not defined, skip site if (args->min_af_type == ALLELE_NONREF && args->min_af>non_ref_ac/(double)an) return 0; // min AF else if (args->min_af_type == ALLELE_MINOR && args->min_af>minor_ac/(double)an) return 0; // min minor AF else if (args->min_af_type == ALLELE_ALT1 && args->min_af>args->ac[1]/(double)an) return 0; // min 1st alternate AF else if (args->min_af_type == ALLELE_MAJOR && args->min_af > major_ac/(double)an) return 0; // min major AF else if (args->min_af_type == ALLELE_NONMAJOR && args->min_af > (an-major_ac)/(double)an) return 0; // min non-major AF } if (args->max_af) { if (an == 0) return 0; // freq not defined, skip site if (args->max_af_type == ALLELE_NONREF && args->max_af<non_ref_ac/(double)an) return 0; // max AF else if (args->max_af_type == ALLELE_MINOR && args->max_af<minor_ac/(double)an) return 0; // max minor AF else if (args->max_af_type == ALLELE_ALT1 && args->max_af<args->ac[1]/(double)an) return 0; // max 1st alternate AF else if (args->max_af_type == ALLELE_MAJOR && args->max_af < major_ac/(double)an) return 0; // max major AF else if (args->max_af_type == ALLELE_NONMAJOR && args->max_af < (an-major_ac)/(double)an) return 0; // max non-major AF } if (args->uncalled) { if (args->uncalled == FLT_INCLUDE && an > 0) return 0; // select uncalled if (args->uncalled == FLT_EXCLUDE && an == 0) return 0; // skip if uncalled } if (args->calc_ac && args->update_info) { bcf_update_info_int32(args->hdr, line, "AC", &args->ac[1], line->n_allele-1); bcf_update_info_int32(args->hdr, line, "AN", &an, 1); } if (args->trim_alts) { int ret = bcf_trim_alleles(args->hsub ? args->hsub : args->hdr, line); if ( ret==-1 ) error("Error: some GT index is out of bounds at %s:%d\n", bcf_seqname(args->hsub ? args->hsub : args->hdr, line), line->pos+1); } if (args->phased) { int phased = bcf_all_phased(args->hdr, line); if (args->phased == FLT_INCLUDE && !phased) { return 0; } // skip unphased if (args->phased == FLT_EXCLUDE && phased) { return 0; } // skip phased } if (args->sites_only) bcf_subset(args->hsub ? args->hsub : args->hdr, line, 0, 0); return 1; }