Back to index

python-biopython  1.60
Static Public Attributes | Private Member Functions
Bio.SeqIO.InsdcIO._InsdcWriter Class Reference
Inheritance diagram for Bio.SeqIO.InsdcIO._InsdcWriter:
Inheritance graph
[legend]
Collaboration diagram for Bio.SeqIO.InsdcIO._InsdcWriter:
Collaboration graph
[legend]

List of all members.

Static Public Attributes

int MAX_WIDTH = 80
int QUALIFIER_INDENT = 21
string QUALIFIER_INDENT_STR = " "
string QUALIFIER_INDENT_TMP = " %s "

Private Member Functions

def _write_feature_qualifier
def _wrap_location
def _write_feature
def _get_annotation_str
def _split_multi_line
def _split_contig

Detailed Description

Base class for GenBank and EMBL writers (PRIVATE).

Definition at line 234 of file InsdcIO.py.


Member Function Documentation

def Bio.SeqIO.InsdcIO._InsdcWriter._get_annotation_str (   self,
  record,
  key,
  default = ".",
  just_first = False 
) [private]
Get an annotation dictionary entry (as a string).

Some entries are lists, in which case if just_first=True the first entry
is returned.  If just_first=False (default) this verifies there is only
one entry before returning it.

Definition at line 310 of file InsdcIO.py.

00310 
00311     def _get_annotation_str(self, record, key, default=".", just_first=False):
00312         """Get an annotation dictionary entry (as a string).
00313 
00314         Some entries are lists, in which case if just_first=True the first entry
00315         is returned.  If just_first=False (default) this verifies there is only
00316         one entry before returning it."""
00317         try:
00318             answer = record.annotations[key]
00319         except KeyError:
00320             return default
00321         if isinstance(answer, list):
00322             if not just_first : assert len(answer) == 1
00323             return str(answer[0])
00324         else:
00325             return str(answer)

def Bio.SeqIO.InsdcIO._InsdcWriter._split_contig (   self,
  record,
  max_len 
) [private]

Definition at line 354 of file InsdcIO.py.

00354 
00355     def _split_contig(self, record, max_len):
00356         "Returns a list of strings, splits on commas."""
00357         #TODO - Merge this with _write_multi_line method?
00358         #It would need the addition of the comma splitting logic...
00359         #are there any other cases where that would be sensible?
00360         contig = record.annotations.get("contig", "")
00361         if isinstance(contig, list) or isinstance(contig, tuple):
00362             contig = "".join(contig)
00363         contig = self.clean(contig)
00364         i = 0
00365         answer = []
00366         while contig:
00367             if len(contig) > max_len:
00368                 #Split lines at the commas
00369                 pos = contig[:max_len-1].rfind(",")
00370                 if pos == -1:
00371                     raise ValueError("Could not break up CONTIG")
00372                 text, contig = contig[:pos+1], contig[pos+1:]
00373             else:
00374                 text, contig = contig, ""
00375             answer.append(text)
00376         return answer

def Bio.SeqIO.InsdcIO._InsdcWriter._split_multi_line (   self,
  text,
  max_len 
) [private]
Returns a list of strings.

Any single words which are too long get returned as a whole line
(e.g. URLs) without an exception or warning.

Definition at line 326 of file InsdcIO.py.

00326 
00327     def _split_multi_line(self, text, max_len):
00328         """Returns a list of strings.
00329         
00330         Any single words which are too long get returned as a whole line
00331         (e.g. URLs) without an exception or warning.
00332         """
00333         #TODO - Do the line spliting while preserving white space?
00334         text = text.strip()
00335         if len(text) <= max_len:
00336             return [text]
00337 
00338         words = text.split()
00339         text = ""
00340         while words and len(text) + 1 + len(words[0]) <= max_len:
00341             text += " " + words.pop(0)
00342             text = text.strip()
00343         #assert len(text) <= max_len
00344         answer = [text]
00345         while words:
00346             text = words.pop(0)
00347             while words and len(text) + 1 + len(words[0]) <= max_len:
00348                 text += " " + words.pop(0)
00349                 text = text.strip()
00350             #assert len(text) <= max_len
00351             answer.append(text)
00352         assert not words
00353         return answer

def Bio.SeqIO.InsdcIO._InsdcWriter._wrap_location (   self,
  location 
) [private]
Split a feature location into lines (break at commas).

Definition at line 275 of file InsdcIO.py.

00275 
00276     def _wrap_location(self, location):
00277         """Split a feature location into lines (break at commas)."""
00278         #TODO - Rewrite this not to recurse!
00279         length = self.MAX_WIDTH - self.QUALIFIER_INDENT
00280         if len(location) <= length:
00281             return location
00282         index = location[:length].rfind(",")
00283         if index == -1:
00284             #No good place to split (!)
00285             import warnings
00286             warnings.warn("Couldn't split location:\n%s" % location)
00287             return location
00288         return location[:index+1] + "\n" + \
00289                self.QUALIFIER_INDENT_STR + self._wrap_location(location[index+1:])

Here is the call graph for this function:

Here is the caller graph for this function:

def Bio.SeqIO.InsdcIO._InsdcWriter._write_feature (   self,
  feature,
  record_length 
) [private]
Write a single SeqFeature object to features table.

Definition at line 290 of file InsdcIO.py.

00290 
00291     def _write_feature(self, feature, record_length):
00292         """Write a single SeqFeature object to features table."""
00293         assert feature.type, feature
00294         location = _insdc_feature_location_string(feature, record_length)
00295         f_type = feature.type.replace(" ","_")
00296         line = (self.QUALIFIER_INDENT_TMP  % f_type)[:self.QUALIFIER_INDENT] \
00297                + self._wrap_location(location) + "\n"
00298         self.handle.write(line)
00299         #Now the qualifiers...
00300         for key, values in feature.qualifiers.iteritems():
00301             if isinstance(values, list) or isinstance(values, tuple):
00302                 for value in values:
00303                     self._write_feature_qualifier(key, value)
00304             elif values:
00305                 #String, int, etc
00306                 self._write_feature_qualifier(key, values)
00307             else:
00308                 #e.g. a /psuedo entry
00309                 self._write_feature_qualifier(key)

Here is the call graph for this function:

def Bio.SeqIO.InsdcIO._InsdcWriter._write_feature_qualifier (   self,
  key,
  value = None,
  quote = None 
) [private]

Definition at line 241 of file InsdcIO.py.

00241 
00242     def _write_feature_qualifier(self, key, value=None, quote=None):
00243         if not value:
00244             self.handle.write("%s/%s\n" % (self.QUALIFIER_INDENT_STR, key))
00245             return
00246         #Quick hack with no line wrapping, may be useful for testing:
00247         #self.handle.write('%s/%s="%s"\n' % (self.QUALIFIER_INDENT_STR, key, value))
00248         if quote is None:
00249             #Try to mimic unwritten rules about when quotes can be left out:
00250             if _is_int_or_long(value):
00251                 quote = False
00252             else:
00253                 quote = True
00254         if quote:
00255             line = '%s/%s="%s"' % (self.QUALIFIER_INDENT_STR, key, value)
00256         else:
00257             line = '%s/%s=%s' % (self.QUALIFIER_INDENT_STR, key, value)
00258         if len(line) <= self.MAX_WIDTH:
00259             self.handle.write(line+"\n")
00260             return
00261         while line.lstrip():
00262             if len(line) <= self.MAX_WIDTH:
00263                 self.handle.write(line+"\n")
00264                 return
00265             #Insert line break...
00266             for index in range(min(len(line)-1, self.MAX_WIDTH),
00267                                self.QUALIFIER_INDENT+1,-1):
00268                 if line[index] == " " : break
00269             if line[index] != " ":
00270                 #No nice place to break...
00271                 index = self.MAX_WIDTH
00272             assert index <= self.MAX_WIDTH
00273             self.handle.write(line[:index] + "\n")
00274             line = self.QUALIFIER_INDENT_STR + line[index:].lstrip()

Here is the call graph for this function:

Here is the caller graph for this function:


Member Data Documentation

Definition at line 236 of file InsdcIO.py.

Reimplemented in Bio.SeqIO.InsdcIO.ImgtWriter, and Bio.SeqIO.InsdcIO.EmblWriter.

Definition at line 238 of file InsdcIO.py.

Reimplemented in Bio.SeqIO.InsdcIO.ImgtWriter, and Bio.SeqIO.InsdcIO.EmblWriter.

Definition at line 239 of file InsdcIO.py.


The documentation for this class was generated from the following file: