2 # Code to read, write and edit VCF files
4 # VCF lines are encoded as a dictionary with these keys (note: all lowercase):
9 # 'alt': list of strings
11 # 'filter': None (missing value), or list of keys (strings); empty list parsed as ["PASS"]
12 # 'info': dictionary of values (see below)
13 # 'format': list of keys (strings)
14 # sample keys: dictionary of values (see below)
16 # The sample keys are accessible through vcf.getsamples()
18 # A dictionary of values contains value keys (defined in ##INFO or ##FORMAT lines) which map
19 # to a list, containign integers, floats, strings, or characters. Missing values are replaced
20 # by a particular value, often -1 or .
22 # Genotypes are not stored as a string, but as a list of 1 or 3 elements (for haploid and diploid samples),
23 # the first (and last) the integer representing an allele, and the second the separation character.
24 # Note that there is just one genotype per sample, but for consistency the single element is stored in a list.
26 # Header lines other than ##INFO, ##FORMAT and ##FILTER are stored as (key, value) pairs and are accessible
29 # The VCF class can be instantiated with a 'regions' variable consisting of tuples (chrom,start,end) encoding
30 # 0-based half-open segments. Only variants with a position inside the segment will be parsed. A regions
31 # parser is available under parse_regions.
33 # When instantiated, a reference can be passed to the VCF class. This may be any class that supports a
34 # fetch(chrom, start, end) method.
38 # NOTE: the position that is returned to Python is 0-based, NOT 1-based as in the VCF file.
43 # only v4.0 writing is complete; alleles are not converted to v3.3 format
46 from collections import namedtuple, defaultdict
47 from operator import itemgetter
48 import sys, re, copy, bisect
54 gtsRegEx = re.compile("[|/\\\\]")
55 alleleRegEx = re.compile('^[ACGTN]+$')
57 # Utility function. Uses 0-based coordinates
58 def get_sequence(chrom, start, end, fa):
59 # obtain sequence from .fa file, without truncation
60 if end<=start: return ""
61 if not fa: return "N"*(end-start)
62 if start<0: return "N"*(-start) + get_sequence(chrom, 0, end, fa).upper()
63 sequence = fa.fetch(chrom, start, end).upper()
64 if len(sequence) < end-start: sequence += "N"*(end-start-len(sequence))
67 # Utility function. Parses a region string
68 def parse_regions( string ):
70 for r in string.split(','):
72 chrom, start, end = elts[0], 0, 3000000000
76 ielts = elts[1].split('-')
77 if len(ielts) != 2: ValueError("Don't understand region string '%s'" % r)
78 try: start, end = int(ielts[0])-1, int(ielts[1])
79 except: raise ValueError("Don't understand region string '%s'" % r)
81 raise ValueError("Don't understand region string '%s'" % r)
82 result.append( (chrom,start,end) )
86 FORMAT = namedtuple('FORMAT','id numbertype number type description missingvalue')
88 ###########################################################################################################
92 ###########################################################################################################
94 cdef class VCFRecord( TabProxies.TupleProxy):
97 initialized from data and vcf meta
104 def __init__(self, vcf):
106 # if len(data) != len(self.vcf._samples):
107 # self.error(str(data),
108 # self.BAD_NUMBER_OF_COLUMNS,
109 # "expected %s for %s samples (%s), got %s" % \
110 # (len(self.vcf._samples),
111 # len(self.vcf._samples),
115 def __cinit__(self, vcf ):
116 # start indexed access at genotypes
121 cdef update( self, char * buffer, size_t nbytes ):
122 '''update internal data.
124 nbytes does not include the terminal '\0'.
126 TabProxies.TupleProxy.update( self, buffer, nbytes )
128 self.contig = self.fields[0]
129 # vcf counts from 1 - correct here
130 self.pos = atoi( self.fields[1] ) - 1
133 return max(0, self.nfields - 9)
136 def __get__(self): return self.contig
139 def __get__(self): return self.pos
142 def __get__(self): return self.fields[2]
145 def __get__(self): return self.fields[3]
149 # convert v3.3 to v4.0 alleles below
151 if alt == ".": alt = []
152 else: alt = alt.upper().split(',')
157 qual = self.fields[5]
158 if qual == b".": qual = -1
160 try: qual = float(qual)
161 except: self.vcf.error(str(self),self.QUAL_NOT_NUMERICAL)
166 # postpone checking that filters exist. Encode missing filter or no filtering as empty list
167 if f == b"." or f == b"PASS" or f == b"0": return []
168 else: return f.split(';')
173 # dictionary of keys, and list of values
176 for blurp in col.split(';'):
177 elts = blurp.split('=')
178 if len(elts) == 1: v = None
179 elif len(elts) == 2: v = elts[1]
180 else: self.vcf.error(str(self),self.ERROR_INFO_STRING)
181 info[elts[0]] = self.vcf.parse_formatdata(elts[0], v, self.vcf._info, str(self))
186 return self.fields[8].split(':')
190 return self.vcf._samples
192 def __getitem__(self, key):
194 # parse sample columns
195 values = self.fields[self.vcf._sample2column[key]].split(':')
199 if len(values) > len(format):
200 self.error(str(self.line),self.BAD_NUMBER_OF_VALUES,"(found %s values in element %s; expected %s)" %\
201 (len(values),key,len(format)))
204 for idx in range(len(format)):
205 expected = self.vcf.get_expected(format[idx], self.vcf._format, alt)
206 if idx < len(values): value = values[idx]
208 if expected == -1: value = "."
209 else: value = ",".join(["."]*expected)
211 result[format[idx]] = self.vcf.parse_formatdata(format[idx], value, self.vcf._format, str(self.data))
212 if expected != -1 and len(result[format[idx]]) != expected:
213 self.error(str(self.data),self.BAD_NUMBER_OF_PARAMETERS,
214 "id=%s, expected %s parameters, got %s" % (format[idx],expected,result[format[idx]]))
215 if len(result[format[idx]] ) < expected: result[format[idx]] += [result[format[idx]][-1]]*(expected-len(result[format[idx]]))
216 result[format[idx]] = result[format[idx]][:expected]
221 cdef class asVCFRecord( ctabix.Parser ):
222 '''converts a :term:`tabix row` into a VCF record.'''
224 def __init__(self, vcffile ):
225 self.vcffile = vcffile
226 def __call__(self, char * buffer, int len ):
228 r = VCFRecord( self.vcffile )
229 r.copy( buffer, len )
240 NT_PHASED_GENOTYPES = 5
242 _errors = { 0:"UNKNOWN_FORMAT_STRING:Unknown file format identifier",
243 1:"BADLY_FORMATTED_FORMAT_STRING:Formatting error in the format string",
244 2:"BADLY_FORMATTED_HEADING:Did not find 9 required headings (CHROM, POS, ..., FORMAT) %s",
245 3:"BAD_NUMBER_OF_COLUMNS:Wrong number of columns found (%s)",
246 4:"POS_NOT_NUMERICAL:Position column is not numerical",
247 5:"UNKNOWN_CHAR_IN_REF:Unknown character in reference field",
248 6:"V33_BAD_REF:Reference should be single-character in v3.3 VCF",
249 7:"V33_BAD_ALLELE:Cannot interpret allele for v3.3 VCF",
250 8:"POS_NOT_POSITIVE:Position field must be >0",
251 9:"QUAL_NOT_NUMERICAL:Quality field must be numerical, or '.'",
252 10:"ERROR_INFO_STRING:Error while parsing info field",
253 11:"ERROR_UNKNOWN_KEY:Unknown key (%s) found in formatted field (info; format; or filter)",
254 12:"ERROR_FORMAT_NOT_NUMERICAL:Expected integer or float in formatted field; got %s",
255 13:"ERROR_FORMAT_NOT_CHAR:Eexpected character in formatted field; got string",
256 14:"FILTER_NOT_DEFINED:Identifier (%s) in filter found which was not defined in header",
257 15:"FORMAT_NOT_DEFINED:Identifier (%s) in format found which was not defined in header",
258 16:"BAD_NUMBER_OF_VALUES:Found too many of values in sample column (%s)",
259 17:"BAD_NUMBER_OF_PARAMETERS:Found unexpected number of parameters (%s)",
260 18:"BAD_GENOTYPE:Cannot parse genotype (%s)",
261 19:"V40_BAD_ALLELE:Bad allele found for v4.0 VCF (%s)",
262 20:"MISSING_REF:Reference allele missing",
263 21:"V33_UNMATCHED_DELETION:Deleted sequence does not match reference (%s)",
264 22:"V40_MISSING_ANGLE_BRACKETS:Format definition is not deliminted by angular brackets",
265 23:"FORMAT_MISSING_QUOTES:Description field in format definition is not surrounded by quotes",
266 24:"V40_FORMAT_MUST_HAVE_NAMED_FIELDS:Fields in v4.0 VCF format definition must have named fields",
267 25:"HEADING_NOT_SEPARATED_BY_TABS:Heading line appears separated by spaces, not tabs",
268 26:"WRONG_REF:Wrong reference %s",
269 27:"ERROR_TRAILING_DATA:Numerical field ('%s') has semicolon-separated trailing data",
270 28:"BAD_CHR_TAG:Error calculating chr tag for %s",
271 29:"ZERO_LENGTH_ALLELE:Found zero-length allele",
272 30:"MISSING_INDEL_ALLELE_REF_BASE:Indel alleles must begin with single reference base"
275 # tag-value pairs; tags are not unique; does not include fileformat, INFO, FILTER or FORMAT fields
278 # version number; 33=v3.3; 40=v4.0
281 # info, filter and format data
286 # header; and required columns
287 _required = ["CHROM","POS","ID","REF","ALT","QUAL","FILTER","INFO","FORMAT"]
291 _ignored_errors = set([11]) # ERROR_UNKNOWN_KEY
292 _warn_errors = set([])
298 # regions to include; None includes everything
306 def __init__(self, _copy=None, reference=None, regions=None, lines=None, leftalign=False):
307 # make error identifiers accessible by name
308 for id in self._errors.keys(): self.__dict__[self._errors[id].split(':')[0]] = id
310 self._leftalign = _copy._leftalign
311 self._header = _copy._header[:]
312 self._version = _copy._version
313 self._info = copy.deepcopy(_copy._info)
314 self._filter = copy.deepcopy(_copy._filter)
315 self._format = copy.deepcopy(_copy._format)
316 self._samples = _copy._samples[:]
317 self._sample2column = copy.deepcopy(_copy._sample2column)
318 self._ignored_errors = copy.deepcopy(_copy._ignored_errors)
319 self._warn_errors = copy.deepcopy(_copy._warn_errors)
320 self._reference = _copy._reference
321 self._regions = _copy._regions
322 if reference: self._reference = reference
323 if regions: self._regions = regions
324 if leftalign: self._leftalign = leftalign
327 def error(self,line,error,opt=None):
328 if error in self._ignored_errors: return
329 errorlabel, errorstring = self._errors[error].split(':')
330 if opt: errorstring = errorstring % opt
331 errwarn = ["Error","Warning"][error in self._warn_errors]
332 sys.stderr.write("Line %s: '%s'\n%s %s: %s\n" % (self._lineno,line,errwarn,errorlabel,errorstring))
333 if error in self._warn_errors: return
334 raise ValueError(errorstring)
336 def parse_format(self,line,format,filter=False):
337 if self._version == 40:
338 if not format.startswith('<'):
339 self.error(line,self.V40_MISSING_ANGLE_BRACKETS)
341 if not format.endswith('>'):
342 self.error(line,self.V40_MISSING_ANGLE_BRACKETS)
344 format = format[1:-1]
345 data = {'id':None,'number':None,'type':None,'descr':None}
347 while len(format.strip())>0:
348 elts = format.strip().split(',')
349 first, rest = elts[0], ','.join(elts[1:])
350 if first.find('=') == -1 or (first.find('"')>=0 and first.find('=') > first.find('"')):
351 if self._version == 40: self.error(line,self.V40_FORMAT_MUST_HAVE_NAMED_FIELDS)
352 if idx == 4: self.error(line,self.BADLY_FORMATTED_FORMAT_STRING)
353 first = ["ID=","Number=","Type=","Description="][idx] + first
354 if first.startswith('ID='): data['id'] = first.split('=')[1]
355 elif first.startswith('Number='): data['number'] = first.split('=')[1]
356 elif first.startswith('Type='): data['type'] = first.split('=')[1]
357 elif first.startswith('Description='):
358 elts = format.split('"')
360 self.error(line,self.FORMAT_MISSING_QUOTES)
361 elts = first.split('=') + [rest]
362 data['descr'] = elts[1]
363 rest = '"'.join(elts[2:])
364 if rest.startswith(','): rest = rest[1:]
366 self.error(line,self.BADLY_FORMATTED_FORMAT_STRING)
369 if filter and idx==1: idx=3 # skip number and type fields for FILTER format strings
370 if not data['id']: self.error(line,self.BADLY_FORMATTED_FORMAT_STRING)
371 if not data['descr']:
372 self.error(line,self.BADLY_FORMATTED_FORMAT_STRING)
373 data['descr'] = '<none>'
374 if not data['type'] and not data['number']:
375 # fine, ##filter format
376 return FORMAT(data['id'],self.NT_NUMBER,0,"Flag",data['descr'],'.')
377 if not data['type'] in ["Integer","Float","Character","String","Flag"]:
378 self.error(line,self.BADLY_FORMATTED_FORMAT_STRING)
379 # I would like a missing-value field, but it isn't there
380 if data['type'] in ['Integer','Float']: data['missing'] = None # Do NOT use arbitrary int/float as missing value
381 else: data['missing'] = '.'
382 if not data['number']: self.error(line,self.BADLY_FORMATTED_FORMAT_STRING)
384 n = int(data['number'])
388 if data['number'] == '.': t = self.NT_UNKNOWN
389 elif data['number'] == '#alleles': t = self.NT_ALLELES
390 elif data['number'] == '#nonref_alleles': t = self.NT_NR_ALLELES
391 elif data['number'] == '#genotypes': t = self.NT_GENOTYPES
392 elif data['number'] == '#phased_genotypes': t = self.NT_PHASED_GENOTYPES
394 self.error(line,self.BADLY_FORMATTED_FORMAT_STRING)
395 return FORMAT(data['id'],t,n,data['type'],data['descr'],data['missing'])
398 def format_format( self, fmt, filter=False ):
399 values = [('ID',fmt.id)]
400 if fmt.number != None and not filter:
401 if fmt.numbertype == self.NT_UNKNOWN: nmb = "."
402 elif fmt.numbertype == self.NT_NUMBER: nmb = str(fmt.number)
403 elif fmt.numbertype == self.NT_ALLELES: nmb = "#alleles"
404 elif fmt.numbertype == self.NT_NR_ALLELES: nmb = "#nonref_alleles"
405 elif fmt.numbertype == self.NT_GENOTYPES: nmb = "#genotypes"
406 elif fmt.numbertype == self.NT_PHASED_GENOTYPES: nmb = "#phased_genotypes"
408 raise ValueError("Unknown number type encountered: %s" % fmt.numbertype)
409 values.append( ('Number',nmb) )
410 values.append( ('Type', fmt.type) )
411 values.append( ('Description', '"' + fmt.description + '"') )
412 if self._version == 33:
413 format = ",".join([v for k,v in values])
415 format = "<" + (",".join( ["%s=%s" % (k,v) for (k,v) in values] )) + ">"
418 def get_expected(self, format, formatdict, alt):
419 fmt = formatdict[format]
420 if fmt.numbertype == self.NT_UNKNOWN: return -1
421 if fmt.numbertype == self.NT_NUMBER: return fmt.number
422 if fmt.numbertype == self.NT_ALLELES: return len(alt)+1
423 if fmt.numbertype == self.NT_NR_ALLELES: return len(alt)
424 if fmt.numbertype == self.NT_GENOTYPES: return ((len(alt)+1)*(len(alt)+2)) // 2
425 if fmt.numbertype == self.NT_PHASED_GENOTYPES: return (len(alt)+1)*(len(alt)+1)
429 def _add_definition(self, formatdict, key, data, line ):
430 if key in formatdict: return
431 self.error(line,self.ERROR_UNKNOWN_KEY,key)
433 formatdict[key] = FORMAT(key,self.NT_NUMBER,0,"Flag","(Undefined tag)",".")
435 if data == []: data = [""] # unsure what type -- say string
436 if type(data[0]) == type(0.0):
437 formatdict[key] = FORMAT(key,self.NT_UNKNOWN,-1,"Float","(Undefined tag)",None)
439 if type(data[0]) == type(0):
440 formatdict[key] = FORMAT(key,self.NT_UNKNOWN,-1,"Integer","(Undefined tag)",None)
442 formatdict[key] = FORMAT(key,self.NT_UNKNOWN,-1,"String","(Undefined tag)",".")
445 # todo: trim trailing missing values
446 def format_formatdata( self, data, format, key=True, value=True, separator=":" ):
447 output, sdata = [], []
448 if type(data) == type([]): # for FORMAT field, make data with dummy values
450 for k in data: d[k] = []
452 # convert missing values; and silently add definitions if required
454 self._add_definition( format, k, data[k], "(output)" )
455 for idx,v in enumerate(data[k]):
456 if v == format[k].missingvalue: data[k][idx] = "."
457 # make sure GT comes first; and ensure fixed ordering; also convert GT data back to string
459 if k != 'GT': sdata.append( (k,data[k]) )
462 sdata = [('GT',map(self.convertGTback,data['GT']))] + sdata
466 if v != None: output.append( k+"="+','.join(map(str,v)) )
467 else: output.append( k )
468 elif key: output.append(k)
470 if v != None: output.append( ','.join(map(str,v)) )
471 else: output.append( "." ) # should not happen
472 # snip off trailing missing data
473 while len(output) > 1:
474 last = output[-1].replace(',','').replace('.','')
475 if len(last)>0: break
477 return separator.join(output)
480 def enter_default_format(self):
481 for f in [FORMAT('GT',self.NT_NUMBER,1,'String','Genotype','.'),
482 FORMAT('GQ',self.NT_NUMBER,1,'Integer','Genotype Quality',-1),
483 FORMAT('DP',self.NT_NUMBER,1,'Integer','Read depth at this position for this sample',-1),
484 FORMAT('HQ',self.NT_UNKNOWN,-1,'Integer','Haplotype Quality',-1), # unknown number, since may be haploid
485 FORMAT('FT',self.NT_NUMBER,1,'String','Sample Genotype Filter','.')]:
486 if f.id not in self._format:
487 self._format[f.id] = f
489 def parse_header( self, line ):
490 assert line.startswith('##')
491 elts = line[2:].split('=')
492 key = elts[0].strip()
493 value = '='.join(elts[1:]).strip()
494 if key == "fileformat":
495 if value == "VCFv3.3":
497 elif value == "VCFv4.0":
499 elif value == "VCFv4.1":
503 self.error(line,self.UNKNOWN_FORMAT_STRING)
505 f = self.parse_format(line, value)
506 self._info[ f.id ] = f
507 elif key == "FILTER":
508 f = self.parse_format(line, value, filter=True)
509 self._filter[ f.id ] = f
510 elif key == "FORMAT":
511 f = self.parse_format(line, value)
512 self._format[ f.id ] = f
514 # keep other keys in the header field
515 self._header.append( (key,value) )
518 def write_header( self, stream ):
519 stream.write("##fileformat=VCFv%s.%s\n" % (self._version // 10, self._version % 10))
520 for key,value in self._header: stream.write("##%s=%s\n" % (key,value))
521 for var,label in [(self._info,"INFO"),(self._filter,"FILTER"),(self._format,"FORMAT")]:
522 for f in var.itervalues(): stream.write("##%s=%s\n" % (label,self.format_format(f,filter=(label=="FILTER"))))
525 def parse_heading( self, line ):
526 assert line.startswith('#')
527 assert not line.startswith('##')
528 headings = line[1:].split('\t')
529 if len(headings)==1 and len(line[1:].split()) >= 9:
530 self.error(line,self.HEADING_NOT_SEPARATED_BY_TABS)
531 headings = line[1:].split()
533 for i,s in enumerate(self._required):
535 if len(headings)<=i or headings[i] != s:
537 if len(headings) <= i:
538 err = "(%sth entry not found)" % (i+1)
540 err = "(found %s, expected %s)" % (headings[i],s)
542 #self.error(line,self.BADLY_FORMATTED_HEADING,err)
544 # allow FORMAT column to be absent
545 if len(headings) == 8:
546 headings.append("FORMAT")
548 self.error(line,self.BADLY_FORMATTED_HEADING,err)
550 self._samples = headings[9:]
551 self._sample2column = dict( [(y,x+9) for x,y in enumerate( self._samples ) ] )
553 def write_heading( self, stream ):
554 stream.write("#" + "\t".join(self._required + self._samples) + "\n")
556 def convertGT(self, GTstring):
557 if GTstring == ".": return ["."]
559 gts = gtsRegEx.split(GTstring)
560 if len(gts) == 1: return [int(gts[0])]
561 if len(gts) != 2: raise ValueError()
562 if gts[0] == "." and gts[1] == ".": return [gts[0],GTstring[len(gts[0]):-len(gts[1])],gts[1]]
563 return [int(gts[0]),GTstring[len(gts[0]):-len(gts[1])],int(gts[1])]
565 self.error(self._line,self.BAD_GENOTYPE,GTstring)
569 def convertGTback(self, GTdata):
570 return ''.join(map(str,GTdata))
572 def parse_formatdata( self, key, value, formatdict, line ):
573 # To do: check that the right number of values is present
574 f = formatdict.get(key,None)
576 self._add_definition(formatdict, key, value, line )
579 if value is not None: self.error(line,self.ERROR_FLAG_HAS_VALUE)
581 values = value.split(',')
582 # deal with trailing data in some early VCF files
583 if f.type in ["Float","Integer"] and len(values)>0 and values[-1].find(';') > -1:
584 self.error(line,self.ERROR_TRAILING_DATA,values[-1])
585 values[-1] = values[-1].split(';')[0]
586 if f.type == "Integer":
587 for idx,v in enumerate(values):
589 if v == ".": values[idx] = f.missingvalue
590 else: values[idx] = int(v)
592 self.error(line,self.ERROR_FORMAT_NOT_NUMERICAL,values)
593 return [0] * len(values)
595 elif f.type == "String":
597 if f.id == "GT": values = map( self.convertGT, values )
599 elif f.type == "Character":
601 if len(v) != 1: self.error(line,self.ERROR_FORMAT_NOT_CHAR)
603 elif f.type == "Float":
604 for idx,v in enumerate(values):
605 if v == ".": values[idx] = f.missingvalue
606 try: return map(float,values)
608 self.error(line,self.ERROR_FORMAT_NOT_NUMERICAL,values)
609 return [0.0] * len(values)
612 self.error(line,self.ERROR_INFO_STRING)
615 def inregion(self, chrom, pos):
616 if not self._regions: return True
617 for r in self._regions:
618 if r[0] == chrom and r[1] <= pos < r[2]: return True
622 def parse_data( self, line, lineparse=False ):
623 cols = line.split('\t')
624 if len(cols) != len(self._samples)+9:
625 # gracefully deal with absent FORMAT column
626 if len(cols) == 8 and len(self._samples)==0:
630 self.BAD_NUMBER_OF_COLUMNS,
631 "expected %s for %s samples (%s), got %s" % (len(self._samples)+9, len(self._samples), self._samples, len(cols)))
635 # get 0-based position
636 try: pos = int(cols[1])-1
637 except: self.error(line,self.POS_NOT_NUMERICAL)
638 if pos < 0: self.error(line,self.POS_NOT_POSITIVE)
640 # implement filtering
641 if not self.inregion(chrom,pos): return None
643 # end of first-pass parse for sortedVCF
644 if lineparse: return chrom, pos, line
648 ref = cols[3].upper()
650 self.error(line,self.MISSING_REF)
651 if self._version == 33: ref = get_sequence(chrom,pos,pos+1,self._reference)
655 if c not in "ACGTN": self.error(line,self.UNKNOWN_CHAR_IN_REF)
656 if "N" in ref: ref = get_sequence(chrom,pos,pos+len(ref),self._reference)
658 # make sure reference is sane
660 left = max(0,pos-100)
661 faref_leftflank = get_sequence(chrom,left,pos+len(ref),self._reference)
662 faref = faref_leftflank[pos-left:]
663 if faref != ref: self.error(line,self.WRONG_REF,"(reference is %s, VCF says %s)" % (faref,ref))
666 # convert v3.3 to v4.0 alleles below
667 if cols[4] == ".": alt = []
668 else: alt = cols[4].upper().split(',')
670 if cols[5] == ".": qual = -1
672 try: qual = float(cols[5])
673 except: self.error(line,self.QUAL_NOT_NUMERICAL)
675 # postpone checking that filters exist. Encode missing filter or no filtering as empty list
676 if cols[6] == "." or cols[6] == "PASS" or cols[6] == "0": filter = []
677 else: filter = cols[6].split(';')
679 # dictionary of keys, and list of values
682 for blurp in cols[7].split(';'):
683 elts = blurp.split('=')
684 if len(elts) == 1: v = None
685 elif len(elts) == 2: v = elts[1]
686 else: self.error(line,self.ERROR_INFO_STRING)
687 info[elts[0]] = self.parse_formatdata(elts[0], v, self._info, line)
689 # Gracefully deal with absent FORMAT column
690 if cols[8] == "": format = []
691 else: format = cols[8].split(':')
693 # check: all filters are defined
695 if f not in self._filter: self.error(line,self.FILTER_NOT_DEFINED, f)
697 # check: format fields are defined
699 if f not in self._format: self.error(line,self.FORMAT_NOT_DEFINED, f)
701 # convert v3.3 alleles
702 if self._version == 33:
703 if len(ref) != 1: self.error(line,self.V33_BAD_REF)
705 have_deletions = False
707 if len(a) == 1: a = a + ref[1:] # SNP; add trailing reference
708 elif a.startswith('I'): a = ref[0] + a[1:] + ref[1:] # insertion just beyond pos; add first and trailing reference
709 elif a.startswith('D'): # allow D<seq> and D<num>
710 have_deletions = True
712 l = int(a[1:]) # throws ValueError if sequence
713 if len(ref) < l: # add to reference if necessary
714 addns = get_sequence(chrom,pos+len(ref),pos+l,self._reference)
716 for i,na in enumerate(newalts): newalts[i] = na+addns
717 a = ref[l:] # new deletion, deleting pos...pos+l
720 if len(ref) < len(s): # add Ns to reference if necessary
721 addns = get_sequence(chrom,pos+len(ref),pos+len(s),self._reference)
722 if not s.endswith(addns) and addns != 'N'*len(addns):
723 self.error(line,self.V33_UNMATCHED_DELETION,
724 "(deletion is %s, reference is %s)" % (a,get_sequence(chrom,pos,pos+len(s),self._reference)))
726 for i,na in enumerate(newalts): newalts[i] = na+addns
727 a = ref[len(s):] # new deletion, deleting from pos
729 self.error(line,self.V33_BAD_ALLELE)
732 # deletion alleles exist, add dummy 1st reference allele, and account for leading base
735 # Petr Danacek's: we can't have a leading nucleotide at (1-based) position 1
736 addn = get_sequence(chrom,pos+len(ref),pos+len(ref)+1,self._reference)
738 alt = [allele+addn for allele in alt]
740 addn = get_sequence(chrom,pos-1,pos,self._reference)
742 alt = [addn + allele for allele in alt]
745 # format v4.0 -- just check for nucleotides
747 if not alleleRegEx.match(allele):
748 self.error(line,self.V40_BAD_ALLELE,allele)
750 # check for leading nucleotide in indel calls
752 if len(allele) != len(ref):
753 if len(allele) == 0: self.error(line,self.ZERO_LENGTH_ALLELE)
754 if ref[0].upper() != allele[0].upper() and "N" not in (ref[0]+allele[0]).upper():
755 self.error(line,self.MISSING_INDEL_ALLELE_REF_BASE)
757 # trim trailing bases in alleles
758 for i in range(1,min(len(ref),min(map(len,alt)))):
759 if len(set(allele[-1].upper() for allele in alt)) > 1 or ref[-1].upper() != alt[0][-1].upper():
761 ref, alt = ref[:-1], [allele[:-1] for allele in alt]
763 # left-align alleles, if a reference is available
764 if self._leftalign and self._reference:
768 if len(allele) > len(ref):
769 longest, shortest = allele, ref
771 longest, shortest = ref, allele
772 if len(longest) == len(shortest) or longest[:len(shortest)].upper() != shortest.upper():
774 if longest[-1].upper() != longest[len(shortest)-1].upper():
779 alt = [allele[:-1] for allele in alt]
780 if min([len(allele) for allele in alt]) == 0 or len(ref) == 0:
781 ref = faref_leftflank[pos-left-1] + ref
782 alt = [faref_leftflank[pos-left-1] + allele for allele in alt]
785 # parse sample columns
787 for sample in cols[9:]:
789 values = sample.split(':')
790 if len(values) > len(format):
791 self.error(line,self.BAD_NUMBER_OF_VALUES,"(found %s values in element %s; expected %s)" % (len(values),sample,len(format)))
792 for idx in range(len(format)):
793 expected = self.get_expected(format[idx], self._format, alt)
794 if idx < len(values): value = values[idx]
796 if expected == -1: value = "."
797 else: value = ",".join(["."]*expected)
798 dict[format[idx]] = self.parse_formatdata(format[idx], value, self._format, line)
799 if expected != -1 and len(dict[format[idx]]) != expected:
800 self.error(line,self.BAD_NUMBER_OF_PARAMETERS,
801 "id=%s, expected %s parameters, got %s" % (format[idx],expected,dict[format[idx]]))
802 if len(dict[format[idx]] ) < expected: dict[format[idx]] += [dict[format[idx]][-1]]*(expected-len(dict[format[idx]]))
803 dict[format[idx]] = dict[format[idx]][:expected]
804 samples.append( dict )
808 'pos':pos, # return 0-based position
816 for key,value in zip(self._samples,samples):
822 def write_data(self, stream, data):
823 required = ['chrom','pos','id','ref','alt','qual','filter','info','format'] + self._samples
825 if k not in data: raise ValueError("Required key %s not found in data" % str(k))
826 if data['alt'] == []: alt = "."
827 else: alt = ",".join(data['alt'])
828 if data['filter'] == None: filter = "."
829 elif data['filter'] == []:
830 if self._version == 33: filter = "0"
831 else: filter = "PASS"
832 else: filter = ';'.join(data['filter'])
833 if data['qual'] == -1: qual = "."
834 else: qual = str(data['qual'])
836 output = [data['chrom'],
837 str(data['pos']+1), # change to 1-based position
843 self.format_formatdata( data['info'], self._info, separator=";" ),
844 self.format_formatdata( data['format'], self._format, value=False ) ]
846 for s in self._samples:
847 output.append( self.format_formatdata( data[s], self._format, key=False ) )
849 stream.write( "\t".join(output) + "\n" )
851 def _parse_header(self, stream):
855 if line.startswith('##'):
856 self.parse_header( line.strip() )
857 elif line.startswith('#'):
858 self.parse_heading( line.strip() )
859 self.enter_default_format()
864 def _parse(self, line, stream):
865 if len(line.strip()) > 0:
866 d = self.parse_data( line.strip() )
870 if self._lines and self._lineno > self._lines: raise StopIteration
871 d = self.parse_data( line.strip() )
874 ######################################################################################################
878 ######################################################################################################
880 def getsamples(self):
881 """ List of samples in VCF file """
884 def setsamples(self,samples):
885 """ List of samples in VCF file """
886 self._samples = samples
889 """ List of header key-value pairs (strings) """
892 def setheader(self,header):
893 """ List of header key-value pairs (strings) """
894 self._header = header
897 """ Dictionary of ##INFO tags, as VCF.FORMAT values """
900 def setinfo(self,info):
901 """ Dictionary of ##INFO tags, as VCF.FORMAT values """
905 """ Dictionary of ##FORMAT tags, as VCF.FORMAT values """
908 def setformat(self,format):
909 """ Dictionary of ##FORMAT tags, as VCF.FORMAT values """
910 self._format = format
913 """ Dictionary of ##FILTER tags, as VCF.FORMAT values """
916 def setfilter(self,filter):
917 """ Dictionary of ##FILTER tags, as VCF.FORMAT values """
918 self._filter = filter
920 def setversion(self, version):
921 if version != 33 and version != 40: raise ValueError("Can only handle v3.3 and v4.0 VCF files")
922 self._version = version
924 def setregions(self, regions):
925 self._regions = regions
927 def setreference(self, ref):
928 """ Provide a reference sequence; a Python class supporting a fetch(chromosome, start, end) method, e.g. PySam.FastaFile """
929 self._reference = ref
931 def ignoreerror(self, errorstring):
932 try: self._ignored_errors.add(self.__dict__[errorstring])
933 except KeyError: raise ValueError("Invalid error string: %s" % errorstring)
935 def warnerror(self, errorstring):
936 try: self._warn_errors.add(self.__dict__[errorstring])
937 except KeyError: raise ValueError("Invalid error string: %s" % errorstring)
939 def parse(self, stream):
940 """ Parse a stream of VCF-formatted lines. Initializes class instance and return generator """
941 last_line = self._parse_header(stream)
942 # now return a generator that does the actual work. In this way the pre-processing is done
943 # before the first piece of data is yielded
944 return self._parse(last_line, stream)
946 def write(self, stream, datagenerator):
947 """ Writes a VCF file to a stream, using a data generator (or list) """
948 self.write_header(stream)
949 self.write_heading(stream)
950 for data in datagenerator: self.write_data(stream,data)
952 def writeheader(self, stream):
953 """ Writes a VCF header """
954 self.write_header(stream)
955 self.write_heading(stream)
957 def compare_calls(self, pos1, ref1, alt1, pos2, ref2, alt2):
958 """ Utility function: compares two calls for equality """
959 # a variant should always be assigned to a unique position, one base before
960 # the leftmost position of the alignment gap. If this rule is implemented
961 # correctly, the two positions must be equal for the calls to be identical.
962 if pos1 != pos2: return False
963 # from both calls, trim rightmost bases when identical. Do this safely, i.e.
964 # only when the reference bases are not Ns
965 while len(ref1)>0 and len(alt1)>0 and ref1[-1] == alt1[-1]:
968 while len(ref2)>0 and len(alt2)>0 and ref2[-1] == alt2[-1]:
971 # now, the alternative alleles must be identical
974 ###########################################################################################################
975 ###########################################################################################################
976 ## API functions added by Andreas
977 ###########################################################################################################
979 def connect( self, filename ):
980 '''connect to tabix file.'''
981 self.tabixfile = pysam.Tabixfile( filename )
982 self._parse_header(self.tabixfile.header)
989 """ Parse a stream of VCF-formatted lines. Initializes class instance and return generator """
991 return self.tabixfile.fetch( reference, start, end, region, parser = asVCFRecord( self ) )
993 def validate( self, record ):
994 '''validate vcf record.
996 returns a validated record.
999 raise NotImplementedError( "needs to be checked" )
1001 chrom, pos = record.chrom, record.pos
1006 self.error(str(record),self.MISSING_REF)
1007 if self._version == 33: ref = get_sequence(chrom,pos,pos+1,self._reference)
1011 if c not in "ACGTN": self.error(str(record),self.UNKNOWN_CHAR_IN_REF)
1012 if "N" in ref: ref = get_sequence(chrom,
1017 # make sure reference is sane
1019 left = max(0,self.pos-100)
1020 faref_leftflank = get_sequence(chrom,left,self.pos+len(ref),self._reference)
1021 faref = faref_leftflank[pos-left:]
1022 if faref != ref: self.error(str(record),self.WRONG_REF,"(reference is %s, VCF says %s)" % (faref,ref))
1025 # check: format fields are defined
1026 for f in record.format:
1027 if f not in self._format: self.error(str(record),self.FORMAT_NOT_DEFINED, f)
1029 # check: all filters are defined
1030 for f in record.filter:
1031 if f not in self._filter: self.error(str(record),self.FILTER_NOT_DEFINED, f)
1033 # convert v3.3 alleles
1034 if self._version == 33:
1035 if len(ref) != 1: self.error(str(record),self.V33_BAD_REF)
1037 have_deletions = False
1039 if len(a) == 1: a = a + ref[1:] # SNP; add trailing reference
1040 elif a.startswith('I'): a = ref[0] + a[1:] + ref[1:] # insertion just beyond pos; add first and trailing reference
1041 elif a.startswith('D'): # allow D<seq> and D<num>
1042 have_deletions = True
1044 l = int(a[1:]) # throws ValueError if sequence
1045 if len(ref) < l: # add to reference if necessary
1046 addns = get_sequence(chrom,pos+len(ref),pos+l,self._reference)
1048 for i,na in enumerate(newalts): newalts[i] = na+addns
1049 a = ref[l:] # new deletion, deleting pos...pos+l
1052 if len(ref) < len(s): # add Ns to reference if necessary
1053 addns = get_sequence(chrom,pos+len(ref),pos+len(s),self._reference)
1054 if not s.endswith(addns) and addns != 'N'*len(addns):
1055 self.error(str(record),self.V33_UNMATCHED_DELETION,
1056 "(deletion is %s, reference is %s)" % (a,get_sequence(chrom,pos,pos+len(s),self._reference)))
1058 for i,na in enumerate(newalts): newalts[i] = na+addns
1059 a = ref[len(s):] # new deletion, deleting from pos
1061 self.error(str(record),self.V33_BAD_ALLELE)
1064 # deletion alleles exist, add dummy 1st reference allele, and account for leading base
1067 # Petr Danacek's: we can't have a leading nucleotide at (1-based) position 1
1068 addn = get_sequence(chrom,pos+len(ref),pos+len(ref)+1,self._reference)
1070 alt = [allele+addn for allele in alt]
1072 addn = get_sequence(chrom,pos-1,pos,self._reference)
1074 alt = [addn + allele for allele in alt]
1077 # format v4.0 -- just check for nucleotides
1079 if not alleleRegEx.match(allele):
1080 self.error(str(record),self.V40_BAD_ALLELE,allele)
1083 # check for leading nucleotide in indel calls
1085 if len(allele) != len(ref):
1086 if len(allele) == 0: self.error(str(record),self.ZERO_LENGTH_ALLELE)
1087 if ref[0].upper() != allele[0].upper() and "N" not in (ref[0]+allele[0]).upper():
1088 self.error(str(record),self.MISSING_INDEL_ALLELE_REF_BASE)
1090 # trim trailing bases in alleles
1091 for i in range(1,min(len(ref),min(map(len,alt)))):
1092 if len(set(allele[-1].upper() for allele in alt)) > 1 or ref[-1].upper() != alt[0][-1].upper():
1094 ref, alt = ref[:-1], [allele[:-1] for allele in alt]
1096 # left-align alleles, if a reference is available
1097 if self._leftalign and self._reference:
1101 if len(allele) > len(ref):
1102 longest, shortest = allele, ref
1104 longest, shortest = ref, allele
1105 if len(longest) == len(shortest) or longest[:len(shortest)].upper() != shortest.upper():
1107 if longest[-1].upper() != longest[len(shortest)-1].upper():
1112 alt = [allele[:-1] for allele in alt]
1113 if min([len(allele) for allele in alt]) == 0 or len(ref) == 0:
1114 ref = faref_leftflank[pos-left-1] + ref
1115 alt = [faref_leftflank[pos-left-1] + allele for allele in alt]
1119 "VCF", "VCFRecord", ]