Closed
Description
The following reproduces the segmentation fault:
$ (echo "##fileformat=VCFv4.2"
echo "##contig=<ID=chr15,length=101991189>"
echo -e "#CHROM\tPOS\tID\tREF\tALT\tQUAL\tFILTER\tINFO"
echo -e "chr15\t23065293\t.\tG\t.\t.\t.\t.") | \
bcftools query -f "%VKX"
Segmentation fault (core dumped)
What is going on is that the conversion function in convert.c
for VKX:
static void process_variantkey_hex(convert_t *convert, bcf1_t *line, fmt_t *fmt, int isample, kstring_t *str)
{
uint64_t vk = variantkey(
convert->header->id[BCF_DT_CTG][line->rid].key,
strlen(convert->header->id[BCF_DT_CTG][line->rid].key),
line->pos,
line->d.allele[0],
strlen(line->d.allele[0]),
line->d.allele[1],
strlen(line->d.allele[1]));
ksprintf(str, "%016" PRIx64 "", vk);
}
does not have a safety check and it attempts to access line->d.allele[1]
even when this pointer is NULL
. This is contrast with the conversion function in convert.c
for ALT:
static void process_alt(convert_t *convert, bcf1_t *line, fmt_t *fmt, int isample, kstring_t *str)
{
int i;
if ( line->n_allele==1 )
{
kputc('.', str);
return;
}
if ( fmt->subscript>=0 )
{
if ( line->n_allele > fmt->subscript+1 )
kputs(line->d.allele[fmt->subscript+1], str);
else
kputc('.', str);
return;
}
for (i=1; i<line->n_allele; i++)
{
if ( i>1 ) kputc(',', str);
kputs(line->d.allele[i], str);
}
}
which properly checks whether line->n_allele==1
instead and indeed:
$ (echo "##fileformat=VCFv4.2"
echo "##contig=<ID=chr15,length=101991189>"
echo -e "#CHROM\tPOS\tID\tREF\tALT\tQUAL\tFILTER\tINFO"
echo -e "chr15\t23065293\t.\tG\t.\t.\t.\t.") | \
bcftools query -f "%ALT"
.
Metadata
Metadata
Assignees
Labels
No labels