Back to index

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

List of all members.

Public Member Functions

def __init__
def write_record
def write_header
def write_footer
def write_records
def write_file
def clean

Public Attributes

 wrap
 record2title
 handle

Private Attributes

 _record_written

Detailed Description

Class to write Fasta format files.

Definition at line 77 of file FastaIO.py.


Constructor & Destructor Documentation

def Bio.SeqIO.FastaIO.FastaWriter.__init__ (   self,
  handle,
  wrap = 60,
  record2title = None 
)
Create a Fasta writer.

handle - Handle to an output file, e.g. as returned
 by open(filename, "w")
wrap -   Optional line length used to wrap sequence lines.
 Defaults to wrapping the sequence at 60 characters
 Use zero (or None) for no wrapping, giving a single
 long line for the sequence.
record2title - Optional function to return the text to be
 used for the title line of each record.  By default the
 a combination of the record.id and record.description
 is used.  If the record.description starts with the
 record.id, then just the record.description is used.

You can either use:

myWriter = FastaWriter(open(filename,"w"))
writer.write_file(myRecords)

Or, follow the sequential file writer system, for example:

myWriter = FastaWriter(open(filename,"w"))
writer.write_header() # does nothing for Fasta files
...
Multiple calls to writer.write_record() and/or writer.write_records()
...
writer.write_footer() # does nothing for Fasta files
writer.close()

Definition at line 79 of file FastaIO.py.

00079 
00080     def __init__(self, handle, wrap=60, record2title=None):
00081         """Create a Fasta writer.
00082 
00083         handle - Handle to an output file, e.g. as returned
00084                  by open(filename, "w")
00085         wrap -   Optional line length used to wrap sequence lines.
00086                  Defaults to wrapping the sequence at 60 characters
00087                  Use zero (or None) for no wrapping, giving a single
00088                  long line for the sequence.
00089         record2title - Optional function to return the text to be
00090                  used for the title line of each record.  By default the
00091                  a combination of the record.id and record.description
00092                  is used.  If the record.description starts with the
00093                  record.id, then just the record.description is used.
00094 
00095         You can either use:
00096 
00097         myWriter = FastaWriter(open(filename,"w"))
00098         writer.write_file(myRecords)
00099 
00100         Or, follow the sequential file writer system, for example:
00101 
00102         myWriter = FastaWriter(open(filename,"w"))
00103         writer.write_header() # does nothing for Fasta files
00104         ...
00105         Multiple calls to writer.write_record() and/or writer.write_records()
00106         ...
00107         writer.write_footer() # does nothing for Fasta files
00108         writer.close()
00109         """
00110         SequentialSequenceWriter.__init__(self, handle)
00111         #self.handle = handle
00112         self.wrap = None
00113         if wrap:
00114             if wrap < 1:
00115                 raise ValueError
00116         self.wrap = wrap
00117         self.record2title = record2title


Member Function Documentation

def Bio.SeqIO.Interfaces.SequenceWriter.clean (   self,
  text 
) [inherited]
Use this to avoid getting newlines in the output.

Definition at line 166 of file Interfaces.py.

00166 
00167     def clean(self, text):
00168         """Use this to avoid getting newlines in the output."""
00169         return text.replace("\n", " ").replace("\r", " ").replace("  ", " ")
    

Here is the caller graph for this function:

def Bio.SeqIO.Interfaces.SequentialSequenceWriter.write_file (   self,
  records 
) [inherited]
Use this to write an entire file containing the given records.

records - A list or iterator returning SeqRecord objects

This method can only be called once.  Returns the number of records
written.

Reimplemented from Bio.SeqIO.Interfaces.SequenceWriter.

Definition at line 262 of file Interfaces.py.

00262 
00263     def write_file(self, records):
00264         """Use this to write an entire file containing the given records.
00265 
00266         records - A list or iterator returning SeqRecord objects
00267 
00268         This method can only be called once.  Returns the number of records
00269         written.
00270         """
00271         self.write_header()
00272         count = self.write_records(records)
00273         self.write_footer()
00274         return count

Here is the call graph for this function:

Definition at line 218 of file Interfaces.py.

00218 
00219     def write_footer(self):
00220         assert self._header_written, "You must call write_header() first"
00221         assert self._record_written, "You have not called write_record() or write_records() yet"
00222         assert not self._footer_written, "You have aleady called write_footer()"
00223         self._footer_written = True

Here is the caller graph for this function:

Definition at line 212 of file Interfaces.py.

00212 
00213     def write_header(self):
00214         assert not self._header_written, "You have aleady called write_header()"
00215         assert not self._record_written, "You have aleady called write_record() or write_records()"
00216         assert not self._footer_written, "You have aleady called write_footer()"
00217         self._header_written = True
        

Here is the caller graph for this function:

def Bio.SeqIO.FastaIO.FastaWriter.write_record (   self,
  record 
)
Write a single Fasta record to the file.

Reimplemented from Bio.SeqIO.Interfaces.SequentialSequenceWriter.

Definition at line 118 of file FastaIO.py.

00118 
00119     def write_record(self, record):
00120         """Write a single Fasta record to the file."""
00121         assert self._header_written
00122         assert not self._footer_written
00123         self._record_written = True
00124 
00125         if self.record2title:
00126             title=self.clean(self.record2title(record))
00127         else:
00128             id = self.clean(record.id)
00129             description = self.clean(record.description)
00130 
00131             #if description[:len(id)]==id:
00132             if description and description.split(None,1)[0]==id:
00133                 #The description includes the id at the start
00134                 title = description
00135             elif description:
00136                 title = "%s %s" % (id, description)
00137             else:
00138                 title = id
00139 
00140         assert "\n" not in title
00141         assert "\r" not in title
00142         self.handle.write(">%s\n" % title)
00143 
00144         data = self._get_seq_string(record) #Catches sequence being None
00145 
00146         assert "\n" not in data
00147         assert "\r" not in data
00148 
00149         if self.wrap:
00150             for i in range(0, len(data), self.wrap):
00151                 self.handle.write(data[i:i+self.wrap] + "\n")
00152         else:
00153             self.handle.write(data + "\n")

def Bio.SeqIO.Interfaces.SequentialSequenceWriter.write_records (   self,
  records 
) [inherited]
Write multiple record to the output file.

records - A list or iterator returning SeqRecord objects

Once you have called write_header() you can call write_record()
and/or write_records() as many times as needed.  Then call
write_footer() and close().

Returns the number of records written.

Definition at line 240 of file Interfaces.py.

00240 
00241     def write_records(self, records):
00242         """Write multiple record to the output file.
00243 
00244         records - A list or iterator returning SeqRecord objects
00245 
00246         Once you have called write_header() you can call write_record()
00247         and/or write_records() as many times as needed.  Then call
00248         write_footer() and close().
00249 
00250         Returns the number of records written.
00251         """
00252         #Default implementation:
00253         assert self._header_written, "You must call write_header() first"
00254         assert not self._footer_written, "You have already called write_footer()"
00255         count = 0
00256         for record in records:
00257             self.write_record(record)
00258             count += 1
00259         #Mark as true, even if there where no records
00260         self._record_written = True
00261         return count

Here is the call graph for this function:

Here is the caller graph for this function:


Member Data Documentation

Reimplemented from Bio.SeqIO.Interfaces.SequentialSequenceWriter.

Definition at line 122 of file FastaIO.py.

Reimplemented from Bio.SeqIO.Interfaces.SequenceWriter.

Definition at line 207 of file Interfaces.py.

Definition at line 116 of file FastaIO.py.

Definition at line 111 of file FastaIO.py.


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