Browse Source

Add support for the Gene Construction Kit format.

Provides a reader for the 'gck' format generated by the
Gene Construction Kit software.
master
Damien Goutte-Gattat 2 years ago
parent
commit
e33d85666a
  1. 2
      README.md
  2. 211
      incenp/bio/seqio/GckIO.py
  3. 7
      incenp/bio/seqio/__init__.py

2
README.md

@ -14,6 +14,8 @@ Formats supported
reading and writing supported
* `snapgene` format, used by [SnapGene](https://www.snapgene.com/):
reading support only
* `gck` format, used by [Gene Construction Kit](http://www.textco.com/gene-construction-kit.php):
reading support only
Usage

211
incenp/bio/seqio/GckIO.py

@ -0,0 +1,211 @@
# -*- coding: utf-8 -*-
# Copyright © 2019 Damien Goutte-Gattat
#
# Redistribution and use of this script, with or without modifications,
# are permitted provided that the following conditions are met:
#
# 1. Redistributions of this script must retain the above copyright
# notice, this list of conditions and the following disclaimer.
#
# THIS SOFTWARE IS PROVIDED BY THE AUTHOR ``AS IS'' AND ANY EXPRESS OR
# IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED
# WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
# ARE DISCLAIMED. IN NO EVENT SHALL THE AUTHOR BE LIABLE FOR ANY
# DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL
# DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
# GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
# INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER
# IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR
# OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN
# IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
"""Gene Construction kit (GCK) format support for Biopython
Provide read support for the binary format used by Gene Construction Kit.
"""
from struct import unpack
from Bio import Alphabet
from Bio.Seq import Seq
from Bio.SeqFeature import SeqFeature, FeatureLocation
from Bio.SeqRecord import SeqRecord
def _read(handle, length):
data = handle.read(length)
if len(data) < length:
raise ValueError("Cannot read {} bytes from handle".format(length))
return data
def _read_packet(handle):
length = _read(handle, 4)
length = unpack('>I', length)[0]
data = _read(handle, length)
return (data, length)
def _read_pstring(handle):
"""Read a Pascal string.
A Pascal string is one byte for length followed by the actual string."""
length = _read(handle, 1)
length = unpack('>B', length)[0]
data = _read(handle, length).decode('ASCII')
return data
def _read_p4string(handle):
"""Read a 32-bit Pascal string.
Similar to a Pascal string but length is encoded on 4 bytes."""
length = _read(handle, 4)
length = unpack('>I', length)[0]
data = _read(handle, length).decode('ASCII')
return data
def GckIterator(handle):
# Skip file header
# GCK files start with a 24-bytes header. Bytes 4 and 8 seem to
# always be 12, maybe this could act as a magic cookie. Bytes
# 17-20 and 21-24 contain variable values of unknown meaning.
_read(handle, 24)
# Read the actual sequence data
packet, length = _read_packet(handle)
# The body of the sequence packet starts with a 32-bit integer
# representing the length of the sequence.
seq_length = unpack('>I', packet[:4])[0]
# This length should not be larger than the length of the
# sequence packet.
if seq_length > length - 4:
raise ValueError("Conflicting sequence length values.")
sequence = packet[4:].decode('ASCII')
record = SeqRecord(Seq(sequence, alphabet=Alphabet.generic_dna))
# Skip unknown packet
_read_packet(handle)
# Read features packet
packet, length = _read_packet(handle)
(seq_length, num_features) = unpack('>IH', packet[:6])
# Check that length in the features packet matches the actual
# length of the sequence
if seq_length != len(record):
raise ValueError("Conflicting sequence length values.")
# Each feature is stored in a 92-bytes structure.
if length - 6 != num_features * 92:
raise ValueError("Features packet size inconsistent with number of features.")
for i in range(0, num_features):
offset = 6 + i * 92
feature_data = packet[offset:offset+92]
# There's probably more stuff to unpack in that structure,
# but those values are the only ones I understand.
(start, end, type, strand, has_name, has_comment, version) = \
unpack('>II6xH14xB17xII35xB', feature_data)
if strand == 1: # Reverse strand
strand = -1
elif strand == 2: # Forward strand
strand = 1
elif strand == 3: # Both strands
# Treated the same as a forward strand as BioPython does
# not seem to support dual-stranded features.
strand = 1
location = FeatureLocation(start, end, strand=strand)
# It looks like any value > 0 indicates a CDS...
if type > 0:
type = 'CDS'
else:
type = 'misc_feature'
# Each feature may have a name and a comment, which are then
# stored immediately after the features packet. Names are
# stored as Pascal strings (1 length byte followed by the
# string itself), comments are stored as "32-bit Pascal strings"
# (4 length bytes followed by the string).
qualifiers = {}
if has_name > 0:
name = _read_pstring(handle)
qualifiers['label'] = [name]
if has_comment > 0:
comment = _read_p4string(handle)
qualifiers['note'] = [comment]
# Each feature may exist in several "versions". We keep only
# the most recent version.
if version > 0:
continue
feature = SeqFeature(location, type=type, qualifiers=qualifiers)
record.features.append(feature)
# Read restriction sites packet
# We are not interested in restriction sites, but we must still read
# that packet so that we can skip the names and comments for each
# site, which are stored after that packet in a similar way as for
# the features above.
packet, length = _read_packet(handle)
(seq_length, num_sites) = unpack('>IH', packet[:6])
# Each site is stored in a 88-bytes structure
if length - 6 != num_sites * 88:
raise ValueError("Sites packet size inconsistent with number of sites")
for i in range(0, num_sites):
offset = 6 + i * 88
site_data = packet[offset:offset+88]
(start, end, has_name, has_comment) = unpack('>II24xII48x', site_data)
# Skip names and comments
if has_name:
_read_pstring(handle)
if has_comment:
_read_p4string(handle)
# Skip unknown packet
_read_packet(handle)
# Next in the file are "version packets".
# However they are not properly formatted "packets" as they are not
# preceded by an integer giving their size. Instead we have a
# short integer indicating how many versions are there, and then
# as many 260-bytes block as we have versions.
num_versions = _read(handle, 2)
num_versions = unpack('>H', num_versions)[0]
versions = _read(handle, num_versions * 260)
for i in range(0, num_versions):
offset = i * 260
version_data = versions[offset:offset+260]
# Each version may have a comment, which is then stored
# after all the "version packets".
has_comment = unpack('>I', version_data[-4:])[0]
if has_comment > 0:
_read_p4string(handle)
# Skip unknown fixed-size block
# Whatever this block contains, it is not preceded by any length
# indicator, so I hope its size is indeed constant in all files...
_read(handle, 706)
# Read the construct's name
name = _read_pstring(handle)
if len(name) > 16 or ' ' in name:
# Store that as the record's description
record.description = name
else:
record.name = name
# Circularity byte
# There may be other flags in that block, but their meaning
# is unknown to me.
flags = _read(handle, 17)
circularity = unpack('>16xB', flags)[0]
if circularity > 0:
record.annotations['topology'] = 'circular'
else:
record.annotations['topology'] = 'linear'
yield record

7
incenp/bio/seqio/__init__.py

@ -29,6 +29,7 @@ from Bio import SeqIO
from . import SnapGeneIO
from . import XdnaIO
from . import GckIO
if not 'xdna' in SeqIO._FormatToIterator:
@ -45,3 +46,9 @@ if not 'snapgene' in SeqIO._FormatToIterator:
if not 'snapgene' in SeqIO._BinaryFormats:
SeqIO._BinaryFormats.append('snapgene')
if not 'gck' in SeqIO._FormatToIterator:
SeqIO._FormatToIterator['gck'] = GckIO.GckIterator
if not 'gck' in SeqIO._BinaryFormats:
SeqIO._BinaryFormats.append('gck')
Loading…
Cancel
Save