Skip to content

Commit

Permalink
Use the warnings module rather than printing/writing to stdout and st…
Browse files Browse the repository at this point in the history
…derr. Based on ideas from Eric Talevich (Bug 2754)
  • Loading branch information
peterc committed Apr 23, 2009
1 parent a277384 commit f2b2125
Show file tree
Hide file tree
Showing 8 changed files with 58 additions and 44 deletions.
5 changes: 3 additions & 2 deletions Bio/PDB/Dice.py
Expand Up @@ -4,7 +4,7 @@
# as part of this package.

import re

import warnings
from Bio.PDB.PDBIO import PDBIO


Expand Down Expand Up @@ -41,7 +41,8 @@ def accept_residue(self, residue):
# skip HETATMS
return 0
if icode!=" ":
print "WARNING: Icode at ", residue.get_id()
warnings.warn("WARNING: Icode at %s" % residue.get_id(),
RuntimeError)
if self.start<=resseq<=self.end:
return 1
return 0
Expand Down
4 changes: 2 additions & 2 deletions Bio/PDB/HSExposure.py
Expand Up @@ -4,7 +4,7 @@
# as part of this package.

from math import pi
import sys
import warnings

from Bio.PDB import *
from AbstractPropertyMap import AbstractPropertyMap
Expand Down Expand Up @@ -201,7 +201,7 @@ def pcb_vectors_pymol(self, filename="hs_exp.py"):
@type filename: string
"""
if len(self.ca_cb_list)==0:
sys.stderr.write("Nothing to draw.\n")
warnings.warn("Nothing to draw.", RuntimeWarning)
return
fp=open(filename, "w")
fp.write("from pymol.cgo import *\n")
Expand Down
11 changes: 7 additions & 4 deletions Bio/PDB/MMCIF2Dict.py
Expand Up @@ -76,7 +76,8 @@ def _make_mmcif_dict(self):
token, value=get_token()
# print token, value
if pos!=nr_fields-1:
print "ERROR: broken name-data pair (data missing)!"
warnings.warn("ERROR: broken name-data pair "
"(data missing)!", RuntimeError)
# The last token was not used, so
# don't set token to None! (this means the
# last parsed token goes through the loop again)
Expand All @@ -87,7 +88,8 @@ def _make_mmcif_dict(self):
# print token, value
mmcif_dict[value]=data
if next_token<4:
print "ERROR: broken name-data pair (name-non data pair)!"
warnings.warn("ERROR: broken name-data pair "
"(name-non data pair)!", RuntimeError)
# print token, value
else:
# get next token
Expand All @@ -102,8 +104,9 @@ def _make_mmcif_dict(self):
token=None
else:
# we found some complete garbage
print "ERROR: broken name-data pair (missing name)!"
print token, value
warnings.warn("ERROR: broken name-data pair "
"(missing name)!\n%s %s" % (token, value),
RuntimeError)
mmcif_dict[None].append(value)
# get next token
token=None
Expand Down
22 changes: 12 additions & 10 deletions Bio/PDB/PDBList.py
Expand Up @@ -35,7 +35,8 @@

__doc__="Access the PDB over the internet (for example to download structures)."

import urllib, re, os, sys
import urllib, re, os
from warnings import warn

class PDBList:
"""
Expand Down Expand Up @@ -147,7 +148,7 @@ def get_all_entries(self):
Returns a list of PDB codes in the index file.
"""
entries = []
print "retrieving index file. Takes about 5 MB."
warnings.warn("retrieving index file. Takes about 5 MB.")
url = urllib.urlopen(self.pdb_server+'/pub/pdb/derived_data/index/entries.idx')
# extract four-letter-codes
entries = map(lambda x: x[:4], \
Expand Down Expand Up @@ -232,11 +233,12 @@ def retrieve_pdb_file(self,pdb_code, obsolete=0, compression='.gz',
# check whether the file exists
if not self.overwrite:
if os.path.exists(final_file):
print "file exists, not retrieved",final_file
warnings.warn("file exists, not retrieved %s" % final_file,
RuntimeError)
return final_file

# Retrieve the file
print 'retrieving',url
warnings.warn('retrieving %s' % url)
lines=urllib.urlopen(url).read()
open(filename,'wb').write(lines)
# uncompress the file
Expand All @@ -260,10 +262,10 @@ def update_pdb(self):

for pdb_code in new+modified:
try:
print 'retrieving %s'%(pdb_code)
warnings.warn('retrieving %s' % pdb_code)
self.retrieve_pdb_file(pdb_code)
except:
print 'error %s'%(pdb_code)
warnings.warn('error %s' % pdb_code, RuntimeError)
# you can insert here some more log notes that
# something has gone wrong.

Expand Down Expand Up @@ -311,14 +313,17 @@ def download_obsolete_entries(self,listfile=None):
def get_seqres_file(self,savefile='pdb_seqres.txt'):
"""Retrieves a (big) file containing all the sequences
of PDB entries and writes it to a file."""
print "retrieving sequence file. Takes about 15 MB."
warnings.warn("retrieving sequence file. Takes about 15 MB.")
url = urllib.urlopen(self.pdb_server+'/pub/pdb/derived_data/pdb_seqres.txt')
file = url.readlines()
open(savefile,'w').writelines(file)



if __name__ == '__main__':

import sys

doc = """PDBList.py
(c) Kristian Rother 2003, Contributed to BioPython
Expand Down Expand Up @@ -367,7 +372,4 @@ def get_seqres_file(self,savefile='pdb_seqres.txt'):
elif re.search('^\d...$',sys.argv[1]):
# get single PDB entry
pl.retrieve_pdb_file(sys.argv[1],pdir=pdb_path)




9 changes: 5 additions & 4 deletions Bio/PDB/PDBParser.py
Expand Up @@ -4,8 +4,7 @@
# as part of this package.

# Python stuff
import sys

import warnings
import numpy

# My stuff
Expand Down Expand Up @@ -242,8 +241,10 @@ def _handle_PDB_exception(self, message, line_counter):
message="%s at line %i." % (message, line_counter)
if self.PERMISSIVE:
# just print a warning - some residues/atoms may be missing
print "PDBConstructionException: %s" % message
print "Exception ignored.\nSome atoms or residues may be missing in the data structure."
warnings.warn("PDBConstructionException: %s\n"
"Exception ignored.\n"
"Some atoms or residues may be missing in the data structure."
% message, RuntimeWarning)
else:
# exceptions are fatal - raise again with new message (including line nr)
raise PDBConstructionException(message)
Expand Down
27 changes: 16 additions & 11 deletions Bio/PDB/StructureBuilder.py
Expand Up @@ -3,14 +3,12 @@
# license. Please see the LICENSE file that should have been included
# as part of this package.

import sys


__doc__="""
Consumer class that builds a Structure object. This is used by
the PDBParser and MMCIFparser classes.
"""

import warnings

# My stuff
# SMCRA hierarchy
Expand Down Expand Up @@ -82,8 +80,9 @@ def init_chain(self, chain_id):
if self.model.has_id(chain_id):
self.chain=self.model[chain_id]
if __debug__:
sys.stderr.write("WARNING: Chain %s is discontinuous at line %i.\n"
% (chain_id, self.line_counter))
warnings.warn("WARNING: Chain %s is discontinuous at line %i."
% (chain_id, self.line_counter),
RuntimeWarning)
else:
self.chain=Chain(chain_id)
self.model.add(self.chain)
Expand Down Expand Up @@ -117,8 +116,10 @@ def init_residue(self, resname, field, resseq, icode):
# There already is a residue with the id (field, resseq, icode).
# This only makes sense in the case of a point mutation.
if __debug__:
sys.stderr.write("WARNING: Residue ('%s', %i, '%s') redefined at line %i.\n"
% (field, resseq, icode, self.line_counter))
warnings.warn("WARNING: Residue ('%s', %i, '%s') "
"redefined at line %i."
% (field, resseq, icode, self.line_counter),
RuntimeWarning)
duplicate_residue=self.chain[res_id]
if duplicate_residue.is_disordered()==2:
# The residue in the chain is a DisorderedResidue object.
Expand Down Expand Up @@ -187,8 +188,11 @@ def init_atom(self, name, coord, b_factor, occupancy, altloc, fullname, serial_n
# name of current atom now includes spaces
name=fullname
if __debug__:
sys.stderr.write("WARNING: atom names %s and %s differ only in spaces at line %i.\n"
% (duplicate_fullname, fullname, self.line_counter))
warnings.warn("WARNING: atom names %s and %s differ "
"only in spaces at line %i."
% (duplicate_fullname, fullname,
self.line_counter),
RuntimeWarning)
atom=self.atom=Atom(name, coord, b_factor, occupancy, altloc, fullname, serial_number)
if altloc!=" ":
# The atom is disordered
Expand All @@ -210,8 +214,9 @@ def init_atom(self, name, coord, b_factor, occupancy, altloc, fullname, serial_n
disordered_atom.disordered_add(duplicate_atom)
residue.flag_disordered()
if __debug__:
sys.stderr.write("WARNING: disordered atom found with blank altloc before line %i.\n"
% self.line_counter)
warnings.warn("WARNING: disordered atom found "
"with blank altloc before line %i.\n"
% self.line_counter, RuntimeWarning)
else:
# The residue does not contain this disordered atom
# so we create a new one.
Expand Down
8 changes: 8 additions & 0 deletions Tests/output/test_PDB
Expand Up @@ -2,9 +2,17 @@ test_PDB
PDBConstructionException: Atom N defined twice in residue <Residue ARG het= resseq=2 icode= > at line 19.
Exception ignored.
Some atoms or residues may be missing in the data structure.
WARNING: disordered atom found with blank altloc before line 31.

WARNING: Residue (' ', 4, ' ') redefined at line 41.
PDBConstructionException: Blank altlocs in duplicate residue SER (' ', 4, ' ') at line 41.
Exception ignored.
Some atoms or residues may be missing in the data structure.
WARNING: Residue (' ', 10, ' ') redefined at line 73.
WARNING: Residue (' ', 14, ' ') redefined at line 104.
WARNING: Residue (' ', 16, ' ') redefined at line 133.
WARNING: Residue (' ', 80, ' ') redefined at line 631.
WARNING: Residue (' ', 81, ' ') redefined at line 644.
PDBConstructionException: Atom O defined twice in residue <Residue HOH het=W resseq=67 icode= > at line 820.
Exception ignored.
Some atoms or residues may be missing in the data structure.
Expand Down
16 changes: 5 additions & 11 deletions Tests/test_PDB.py
Expand Up @@ -10,8 +10,10 @@
raise MissingExternalDependencyError(\
"Install NumPy if you want to use Bio.PDB.")

import sys
# Redirect stderr so user does not see warnings
import warnings
def send_warnings_to_stdout(message, category, filename, lineno, file=None):
print message
warnings.showwarning = send_warnings_to_stdout

def quick_neighbor_search_test() :
#Based on the self test in Bio.PDB.NeighborSearch
Expand Down Expand Up @@ -89,12 +91,4 @@ def run_test():
print "NeighborSearch test"
quick_neighbor_search_test()



old_stderr = sys.stderr
# Hide stderr output for user
sys.stderr=TheVoid()
try:
run_test()
finally:
sys.stderr = old_stderr
run_test()

0 comments on commit f2b2125

Please sign in to comment.