FreeSASA  1.1
Open source SASA calculations
View on GitHub
freesasa.pyx
1 ## @package freesasa
2 # @author Simon Mitternacht
3 # @copyright [MIT License](md_license.html)
4 #
5 # @brief FreeSASA Python interface.
6 #
7 # A minimal program would be something like
8 #
9 # ~~~{.py}
10 # structure = freesasa.Structure("1abc.pdb")
11 # result = freesasa.calc(structure)
12 # print "SASA = %f" % result.totalArea()
13 # ~~~
14 #
15 # See documentation of the classes and functions for how to customize behavior
16 #
17 
18 from libc.stdio cimport FILE, fopen, fclose
19 from libc.stdlib cimport free, realloc, malloc
20 from libc.string cimport memcpy
21 from cfreesasa cimport *
22 
23 ## Used to specify the algorithm by Shrake & Rupley
24 ShrakeRupley = 'ShrakeRupley'
25 
26 ## Used to specify the algorithm by Lee & Richards
27 LeeRichards = 'LeeRichards'
28 
29 ## Used for classification
30 polar = 'Polar'
31 
32 ## Used for classification
33 apolar = 'Apolar'
34 
35 ## int: Suppress all warnings and errors (used by setVerbosity())
36 silent = FREESASA_V_SILENT
37 
38 ## int: Suppress all warnings but not errors (used by setVerbosity())
39 nowarnings = FREESASA_V_NOWARNINGS
40 
41 ## int: Normal verbosity (used by setVerbosity())
42 normal = FREESASA_V_NORMAL
43 
44 ## int: Print debug messages (used by setVerbosity())
45 debug = FREESASA_V_DEBUG
46 
47 ## The default values for calculation parameters
48 defaultParameters = {
49  'algorithm' : LeeRichards,
50  'probe-radius' : freesasa_default_parameters.probe_radius,
51  'n-points' : freesasa_default_parameters.shrake_rupley_n_points,
52  'n-slices' : freesasa_default_parameters.lee_richards_n_slices,
53  'n-threads' : freesasa_default_parameters.n_threads
54 }
55 
56 ## Stores parameter values to be used by calculation.
57 #
58 # Wraps the C struct freesasa_parameters.
59 cdef class Parameters:
60 
61  cdef freesasa_parameters _c_param
62  ## Initializes Parameters object.
63  #
64  # @param param (dict) optional argument to specify parameter-values, modeled after
65  # ::defaultParameters
66  # @exception AssertionError invalid parameter values supplied
67  # @see defaultParameters
68  def __init__(self,param=None):
69 
70  self._c_param = freesasa_default_parameters
71  if param != None:
72  if 'algorithm' in param: self.setAlgorithm(param['algorithm'])
73  if 'probe-radius' in param: self.setProbeRadius(param['probe-radius'])
74  if 'n-points' in param: self.setNPoints(param['n-points'])
75  if 'n-slices' in param: self.setNSlices(param['n-slices'])
76  if 'n-threads' in param: self.setNThreads(param['n-threads'])
77  unknownKeys = []
78  for key in param:
79  if not key in defaultParameters:
80  unknownKeys.append(key)
81  if len(unknownKeys) > 0:
82  raise AssertionError('Key(s): ',unknownKeys,', unknown');
83 
84  ## Set algorithm.
85  #
86  # @param alg (str) algorithm name, only allowed values are ::ShrakeRupley and ::LeeRichards
87  # @exception AssertionError unknown algorithm specified
88  def setAlgorithm(self,alg):
89  if alg == ShrakeRupley:
90  self._c_param.alg = FREESASA_SHRAKE_RUPLEY
91  elif alg == LeeRichards:
92  self._c_param.alg = FREESASA_LEE_RICHARDS
93  else:
94  raise AssertionError("Algorithm '%s' is unknown" % alg)
95 
96  ## Get algorithm.
97  #
98  # @return Name of algorithm
99  def algorithm(self):
100  if self._c_param.alg == FREESASA_SHRAKE_RUPLEY:
101  return ShrakeRupley
102  if self._c_param.alg == FREESASA_LEE_RICHARDS:
103  return LeeRichards
104  raise Exception("No algorithm specified, shouldn't be possible")
105 
106  ## Set probe radius.
107  # @param r probe radius in Å (>= 0)
108  # @exception AssertionError r < 0
109  def setProbeRadius(self,r):
110  assert(r >= 0)
111  self._c_param.probe_radius = r
112 
113  ## Get probe radius.
114  # @return Probe radius in Å
115  def probeRadius(self):
116  return self._c_param.probe_radius
117 
118  ## Set number of test points in Shrake & Rupley algorithm.
119  # @param n (int) Number of points (> 0).
120  # @exception AssertionError n <= 0.
121  def setNPoints(self,n):
122  assert(n > 0)
123  self._c_param.shrake_rupley_n_points = n
124 
125  ## Get number of test points in Shrake & Rupley algorithm.
126  # @return Number of points.
127  def nPoints(self):
128  return self._c_param.shrake_rupley_n_points
129 
130  ## Set the number of slices per atom in Lee & Richards algorithm.
131  # @param n (int) Number of slices (> 0)
132  # @exception AssertionError n <= 0
133  def setNSlices(self,n):
134  assert(n> 0)
135  self._c_param.lee_richards_n_slices = n
136 
137  ## Get the number of slices per atom in Lee & Richards algorithm.
138  # @return Number of slices.
139  def nSlices(self):
140  return self._c_param.lee_richards_n_slices
141 
142 
143  ## Set the number of threads to use in calculations.
144  # @param n (int) Number of points (> 0)
145  # @exception AssertionError n <= 0
146  def setNThreads(self,n):
147  assert(n>0)
148  self._c_param.n_threads = n
149 
150 
151  ## Get the number of threads to use in calculations.
152  # @return Number of threads.
153  def nThreads(self):
154  return self._c_param.n_threads
155 
156  # not pretty, but only way I've found to pass pointers around
157  def _get_address(self, size_t ptr2ptr):
158  cdef freesasa_parameters **p = <freesasa_parameters**> ptr2ptr
159  p[0] = &self._c_param
160 
161 ## Stores results from SASA calculation.
162 # The type of object returned by calc(), not intended to be used
163 # outside of that context.
164 #
165 # Wraps the C struct freesasa_result.
166 cdef class Result:
167  cdef freesasa_result* _c_result
168 
169  ## The constructor
170  def __cinit__ (self):
171  self._c_result = NULL
172 
173  ## The destructor
174  def __dealloc__(self):
175  if self._c_result is not NULL:
177 
178  ## Number of atoms in the results
179  # @return Number of atoms
180  def nAtoms(self):
181  if self._c_result is not NULL:
182  return self._c_result.n_atoms
183  return 0
184 
185 
186  ## Total SASA.
187  # @return The total area in Å^2.
188  # @exception AssertionError If no results have been associated
189  # with the object
190  def totalArea(self):
191  assert(self._c_result is not NULL)
192  return self._c_result.total
193 
194  ## SASA for a given atom.
195  # @param i (int) index of atom.
196  # @return SASA of atom i in Å^2.
197  # @exception AssertionError If no results have been associated
198  # with the object or if index is out of bounds
199  def atomArea(self,i):
200  assert(self._c_result is not NULL)
201  assert(i < self._c_result.n_atoms)
202  return self._c_result.sasa[i]
203 
204  def _get_address(self, size_t ptr2ptr):
205  cdef freesasa_result **p = <freesasa_result**> ptr2ptr
206  p[0] = self._c_result
207 
208 
209 ## Assigns class and radius to atom by residue and atom name.
210 #
211 # Subclasses derived from Classifier can be used to define custom
212 # atomic radii and/or classes. Can also be initialized from a @ref
213 # Config-file "configuration file" with a custom classifier.
214 #
215 # Wraps a C freesasa_classifier. If initialized without arguments the
216 # default classifier is used.
217 #
218 # Residue names should be of the format `"ALA"`,`"ARG"`, etc.
219 #
220 # Atom names should be of the format `"CA"`, `"N"`, etc.
221 #
222 cdef class Classifier:
223  cdef freesasa_classifier* _c_classifier
224 
225  ## Constructor.
226  #
227  # If no file is provided the default classifier is used.
228  #
229  # @see @ref Config-file.
230  #
231  # @param fileName Name of file with classifier configuration.
232  # @exception IOError Problem opening/reading file
233  # @exception Exception Problem parsing provided configuration or
234  # initializing defaults
235  def __cinit__ (self,fileName=None):
236  cdef FILE *config
237  self._c_classifier = NULL
238  if fileName is not None:
239  config = fopen(fileName,'r')
240  if config is NULL:
241  raise IOError("File '%s' could not be opened." % fileName)
243  fclose(config)
244  if self._c_classifier is NULL:
245  raise Exception("Error parsing configuration in '%s'." % fileName)
246  else:
247  self._c_classifier = &freesasa_default_classifier
248  ## The destructor
249  def __dealloc__(self):
250  if self._c_classifier is not &freesasa_default_classifier:
252 
253  ## Class of atom.
254  #
255  # Depending on the configuration these classes can be
256  # anything, but typically they will be 'Polar' and 'Apolar'.
257  # Unrecognized atoms will get the class 'Unknown'.
258  #
259  # @param residueName (str) Residue name (`"ALA"`,`"ARG"`,...).
260  # @param atomName (str) Atom name (`"CA"`,`"C"`,...).
261  # @return A string describing the class
262  def classify(self,residueName,atomName):
263  classIndex = self._c_classifier.sasa_class(residueName,atomName,self._c_classifier)
264  if classIndex is FREESASA_WARN:
265  return 'Unknown'
266  return self._c_classifier.class2str(classIndex,self._c_classifier)
267 
268  ## Radius of atom.
269  #
270  # This allows the classifier to be used to calculate the atomic
271  # radii used in calculations. Unknown atoms will get a negative
272  # radius.
273  #
274  # @param residueName (str) Residue name (`"ALA"`, `"ARG"`, ...).
275  # @param atomName (str) Atom name (`"CA"`, `"C"`, ...).
276  # @return The radius in Å.
277  def radius(self,residueName,atomName):
278  return self._c_classifier.radius(residueName,atomName,self._c_classifier)
279 
280  def _get_address(self, size_t ptr2ptr):
281  cdef freesasa_classifier **p = <freesasa_classifier**> ptr2ptr
282  p[0] = self._c_classifier
283 
284 ## Represents a protein structure, including its atomic radii.
285 #
286 # Initialized from PDB-file. Calculates atomic radii using default
287 # classifier, or custom one provided as argument to initalizer
288 #
289 # Wraps the C struct freesasa_structure.
290 #
291 # Since it is intended to be a static structure the word 'get' is
292 # omitted in the getter-functions.
293 cdef class Structure:
294  cdef freesasa_structure* _c_structure
295  ## By default ignore HETATM, Hydrogens, only use first model. For unknown atoms
296  # try to guess the radius, if this fails, assign radius 0 (to
297  # allow changing the radius later).
298  defaultOptions = { 'hetatm' : False,
299  'hydrogen' : False,
300  'join-models' : False,
301  'skip-unknown' : False,
302  'halt-at-unknown' : False
303  }
304 
305  ## Constructor
306  #
307  # If PDB file is provided, the structure will be constructed
308  # based on the file. If not, this simply initializes an empty
309  # structure and the other arguments are ignored. In this case
310  # atoms will have to be added manually using addAtom().
311  #
312  # @param fileName (str) PDB file (if None empty structure generated).
313  # @param classifier An optional Classifier to calculate atomic
314  # radii, uses default if none provided
315  # @param options specify which atoms and models to include
316  # @exception IOError Problem opening/reading file.
317  # @exception Exception Problem parsing PDB file or calculating
318  # atomic radii.
319  # @exception Exception If option 'halt-at-unknown' selected and
320  # unknown atom encountered.
321  def __init__(self,fileName=None,classifier=None,
322  options = defaultOptions):
323 
324  self._c_structure = NULL
325 
326  if fileName is None:
328  return
329  cdef FILE *input
330  input = fopen(fileName,'r')
331  if input is NULL:
332  raise IOError("File '%s' could not be opened." % fileName)
333  structure_options = Structure._get_structure_options(options)
334  self._c_structure = freesasa_structure_from_pdb(input,NULL,structure_options)
335  fclose(input)
336 
337  if self._c_structure is NULL:
338  raise Exception("Error reading '%s'." % fileName)
339 
340  # this means we might be calculating the radii twice, the
341  # advantage of doing it this way is that we can define new
342  # classifiers using a Python interface
343  if (classifier is not None):
344  self.setRadiiWithClassifier(classifier)
345 
346 
347  ## Add atom to structure.
348  #
349  # This function is meant to be used if the structure was not
350  # initialized from a PDB. Default radii will be assigned to each
351  # atom. This can be overriden by calling
352  # setRadiiWithClassifier() afterwards.
353  #
354  # There are no restraints on string lengths for the arguments, but
355  # the atom won't be added if the default classifier doesn't
356  # recognize the atom and also cannot deduce its element from the
357  # atom name.
358  #
359  # @param atomName (str) atom name (e.g. `"CA"`)
360  # @param residueName (str) residue name (e.g. `"ALA"`)
361  # @param residueNumber (str or int) residue number (e.g. `'12'`)
362  # or integer. Some PDBs have residue-numbers that aren't
363  # regular numbers. Therefore treated as a string primarily.
364  # @param chainLabel (str) 1-character string with chain label (e.g. 'A')
365  # @param x,y,z (float) coordinates
366  #
367  # @exception Exception Residue-number invalid
368  def addAtom(self, atomName, residueName, residueNumber, chainLabel, x, y, z):
369  if (type(residueNumber) is str):
370  resnum = residueNumber
371  elif (type(residueNumber) is int):
372  resnum = "%d" % residueNumber
373  else:
374  raise Exception("Residue-number invalid, must be either string or number")
375  cdef const char *label = chainLabel
376  ret = freesasa_structure_add_atom(self._c_structure, atomName,
377  residueName, resnum, label[0],
378  x, y, z)
379  assert(ret != FREESASA_FAIL)
380 
381  ## Assign radii to atoms in structure using a classifier.
382  #
383  # @param classifier A classifier to use to calculate radii
384  # @exception AssertionError if structure not properly initialized
385  def setRadiiWithClassifier(self,classifier):
386  assert(self._c_structure is not NULL)
387  n = self.nAtoms()
388  r = []
389  for i in range(0,n):
390  r.append(classifier.radius(self.residueName(i), self.atomName(i)))
391  self.setRadii(r)
392 
393  ## Set atomic radii from an array
394  # @param radiusArray: Array of atomic radii in Ångström, should
395  # have nAtoms() elements.
396  # @exception AssertionError if radiusArray has wrong dimension, structure
397  # not properly initialized, or if the array contains
398  # negative radii (not properly classified?)
399  def setRadii(self,radiusArray):
400  assert(self._c_structure is not NULL)
401  n = self.nAtoms()
402  assert len(radiusArray) == n
403  cdef double *r = <double *>malloc(sizeof(double)*n)
404  assert(r is not NULL)
405  for i in range(0,n):
406  r[i] = radiusArray[i]
407  assert(r[i] >= 0)
409 
410  ## Number of atoms.
411  #
412  # @return Number of atoms
413  # @exception AssertionError if not properly initialized
414  def nAtoms(self):
415  assert(self._c_structure is not NULL)
417 
418  ## Radius of atom.
419  # @param i (int) Index of atom.
420  # @return Radius in Å.
421  # @exception AssertionError if index out of bounds, object not properly initalized.
422  def radius(self,i):
423  assert(i >= 0 and i < self.nAtoms())
424  assert(self._c_structure is not NULL)
425  cdef const double *r = freesasa_structure_radius(self._c_structure)
426  assert(r is not NULL)
427  return r[i]
428 
429  ## Get atom name
430  # @param i (int) Atom index.
431  # @return Atom name as 4-character string.
432  # @exception AssertionError: if index out of range or Structure not properly initialized.
433  def atomName(self,i):
434  assert(i >= 0 and i < self.nAtoms())
435  assert(self._c_structure is not NULL)
437 
438  ## Get residue name of given atom.
439  # @param i (int) Atom index.
440  # @return Residue name as 3-character string.
441  # @exception AssertionError if index out of range or Structure not properly initialized
442  def residueName(self,i):
443  assert(i >= 0 and i < self.nAtoms())
444  assert(self._c_structure is not NULL)
446 
447  ## Get residue number for given atom
448  # @param i (int) Atom index.
449  # @return Residue number as 4-character string
450  # @exception AssertionError if index out of range or Structure not properly initialized
451  def residueNumber(self,i):
452  assert(i >= 0 and i < self.nAtoms())
453  assert(self._c_structure is not NULL)
455 
456  ## Get chain label for given atom.
457  # @param i (int) Atom index.
458  # @return Chain label as 1-character string
459  # @exception AssertionError if index out of range or Structure not properly initialized
460  def chainLabel(self,i):
461  assert(i >= 0 and i < self.nAtoms())
462  assert(self._c_structure is not NULL)
463  cdef char label[2]
464  label[0] = freesasa_structure_atom_chain(self._c_structure,i);
465  label[1] = '\0'
466  return label
467 
468  ## Get coordinates of given atom.
469  # @param i (int) Atom index.
470  # @return array of x, y, and z coordinates
471  # @exception AssertionError if index out of range or Structure not properly initialized
472  def coord(self, i):
473  assert(i >= 0 and i < self.nAtoms())
474  assert(self._c_structure is not NULL)
475  cdef const double *coord = freesasa_structure_coord_array(self._c_structure);
476  return [coord[3*i], coord[3*i+1], coord[3*i+2]]
477 
478  @staticmethod
479  def _get_structure_options(param):
480  options = 0
481 
482  # check validity of options
483  knownOptions = {'hetatm','hydrogen','join-models','separate-models',
484  'separate-chains','skip-unknown','halt-at-unknown'}
485  unknownOptions = []
486  for key in param:
487  if not key in knownOptions:
488  unknownOptions.append(key)
489  if len(unknownOptions) > 0:
490  raise AssertionError("Option(s): ",unknownOptions," unknown.");
491 
492  # calculate bitfield
493  if 'hetatm' in param and param['hetatm']:
494  options |= FREESASA_INCLUDE_HETATM
495  if 'hydrogen' in param and param['hydrogen']:
496  options |= FREESASA_INCLUDE_HYDROGEN
497  if 'join-models' in param and param['join-models']:
498  options |= FREESASA_JOIN_MODELS
499  if 'separate-models' in param and param['separate-models']:
500  options |= FREESASA_SEPARATE_MODELS
501  if 'separate-chains' in param and param['separate-chains']:
502  options |= FREESASA_SEPARATE_CHAINS
503  if 'skip-unknown' in param and param['skip-unknown']:
504  options |= FREESASA_SKIP_UNKNOWN
505  if 'halt-at-unknown' in param and param['halt-at-unknown']:
506  options |= FREESASA_HALT_AT_UNKNOWN
507  return options
508 
509  def _get_address(self, size_t ptr2ptr):
510  cdef freesasa_structure **p = <freesasa_structure**> ptr2ptr
511  p[0] = self._c_structure
512 
513  def _set_address(self, size_t ptr2ptr):
514  cdef freesasa_structure **p = <freesasa_structure**> ptr2ptr
515  self._c_structure = p[0]
516 
517  ## The destructor
518  def __dealloc__(self):
519  if self._c_structure is not NULL:
521 
522 ## Default options for structureArray.
523 # Defined separately for Doxygen's sake.
524 defaultStructureArrayOptions = {
525  'hetatm' : False,
526  'hydrogen' : False,
527  'separate-chains' : True,
528  'separate-models' : False
529 }
530 
531 ## Create array of structures from PDB file.
532 #
533 # Split PDB file into several structures by either by treating
534 # chains separately, by treating each MODEL as a separate
535 # structure, or both.
536 #
537 # @param fileName (str) The PDB file.
538 # @param options (dict) Specification for how to read the PDB-file
539 # (see def value for options).
540 # @param classifier: Classifier to assign atoms radii, default is used
541 # if none specified.
542 # @exception AssertionError if fileName is None
543 # @exception AssertionError if an option value is not recognized
544 # @exception AssertionError if neither of the options 'separate-chains'
545 # and 'separate-models' are specified.
546 # @exception IOError if can't open file
547 # @exception Exception: if there are problems parsing the input
548 # @return An array of Structures
549 def structureArray(fileName,
550  options = defaultStructureArrayOptions,
551  classifier = None):
552  assert fileName is not None
553  # we need to have at least one of these
554  assert(('separate-chains' in options and options['separate-chains'] is True)
555  or ('separate-models' in options and options['separate-models'] is True))
556  structure_options = Structure._get_structure_options(options)
557  cdef FILE *input
558  input = fopen(fileName,'r')
559  if input is NULL:
560  raise IOError("File '%s' could not be opened." % fileName)
561  cdef int n;
562  cdef freesasa_structure** sArray = freesasa_structure_array(input,&n,NULL,structure_options)
563  fclose(input)
564  if sArray is NULL:
565  raise Exception("Problems reading structures in '%s'." % fileName)
566  structures = []
567  for i in range(0,n):
568  structures.append(Structure())
569  structures[-1]._set_address(<size_t> &sArray[i])
570  if classifier is not None:
571  structures[-1].setRadiiWithClassifier(classifier)
572  free(sArray)
573  return structures
574 
575 
576 ## Calculate SASA of Structure
577 # @param structure Structure to be used
578 # @param parameters Parameters to use (if not specified defaults are used)
579 # @return A Result object
580 # @exception Exception something went wrong in calculation (see C library error messages)
581 def calc(structure,parameters=None):
582  cdef const freesasa_parameters *p = NULL
583  cdef const freesasa_structure *s = NULL
584  if parameters is not None: parameters._get_address(<size_t>&p)
585  structure._get_address(<size_t>&s)
586  result = Result()
587  result._c_result = <freesasa_result*> freesasa_calc_structure(s,p)
588  if result._c_result is NULL:
589  raise Exception("Error calculating SASA.")
590  return result
591 
592 ## Break SASA result down into classes.
593 # @param result Result from sasa calculation.
594 # @param structure Structure used in calculation.
595 # @param classifier Classifier (if not specified default is used).
596 # @return Dictionary with names of classes as keys and their SASA values as values.
597 # @exception Exception: Problems with classification, see C library error messages
598 # (or Python exceptions if run with derived classifier).
599 def classifyResults(result,structure,classifier=None):
600  if classifier is None:
601  classifier = Classifier()
602  ret = dict()
603  for i in range(0,structure.nAtoms()):
604  name = classifier.classify(structure.residueName(i),structure.atomName(i))
605  if name not in ret:
606  ret[name] = 0
607  ret[name] += result.atomArea(i)
608  return ret
609 
610 ## Sum SASA result over a selection of atoms
611 # @param commands A list of commands with selections using Pymol
612 # syntax, e.g. `"s1, resn ala+arg"` or `"s2, chain A and resi 1-5"`
613 # (see @ref Selection).
614 # @param structure A Structure.
615 # @param result Result from sasa calculation on structure.
616 # @return Dictionary with names of selections (`"s1"`,`"s2"`,...) as
617 # keys, and the corresponding SASA values as values.
618 # @exception Exception: Parser failed (typically syntax error), see
619 # library error messages.
620 def selectArea(commands,structure,result):
621  cdef freesasa_structure *s
622  cdef freesasa_result *r
623  cdef double area
624  cdef char *name = <char*>malloc(FREESASA_MAX_SELECTION_NAME+1);
625  structure._get_address(<size_t> &s)
626  result._get_address(<size_t> &r)
627  value = dict()
628  for cmd in commands:
629  ret = freesasa_select_area(cmd,name,&area,s,r)
630  if ret == FREESASA_FAIL:
631  raise Exception("Error parsing '%s'" % cmd)
632  value[name] = area
633  free(name)
634  return value
635 
636 ## Set global verbosity
637 # @param verbosity Can have values freesasa.silent, freesasa.nowarnings or freesasa.normal
638 # @exception AssertionError if verbosity has illegal value
639 def setVerbosity(verbosity):
640  assert(verbosity in [silent, nowarnings, normal])
641  freesasa_set_verbosity(verbosity)
642 
643 
644 ## Get global verbosity
645 # @return Verbosity (freesasa.silent, freesasa.nowarnings or freesasa.normal)
647  return freesasa_get_verbosity()
648 
649 ## Create a freesasa structure from a Bio.PDB structure
650 #
651 # @remark Experimental, not thorougly tested yet
652 #
653 # @param bioPDBStructure a Bio.PDB structure
654 # @param classifier an optional classifier to specify atomic radii
655 # @param options Options supported are 'hetatm', 'skip-unknown' and 'halt-at-unknown'
656 # @return a freesasa.Structure
657 #
658 # @exception Exception if option 'halt-at-unknown' is selected and
659 # unknown atoms are encountered. Passes on exceptions from
660 # Structure.addAtom() and
661 # Structure.setRadiiWithClassifier().
662 def structureFromBioPDB(bioPDBStructure, classifier=None, options = Structure.defaultOptions):
663  structure = Structure()
664  if (classifier is None):
665  classifier = Classifier()
666  optbitfield = Structure._get_structure_options(options)
667 
668  atoms = bioPDBStructure.get_atoms()
669 
670  for a in atoms:
671  r = a.get_parent()
672  hetflag, resseq, icode = r.get_id()
673 
674  if (hetflag is not ' ' and not (optbitfield & FREESASA_INCLUDE_HETATM)):
675  continue
676 
677  c = r.get_parent()
678  v = a.get_vector()
679 
680  if (classifier.classify(r.get_resname(), a.get_fullname()) is 'Unknown'):
681  if (optbitfield & FREESASA_SKIP_UNKNOWN):
682  continue
683  if (optbitfield & FREESASA_HALT_AT_UNKNOWN):
684  raise Exception("Halting at unknown atom")
685 
686  structure.addAtom(a.get_fullname(), r.get_resname(), resseq, c.get_id(),
687  v[0], v[1], v[2])
688 
689  structure.setRadiiWithClassifier(classifier)
690  return structure
691 
692 ## Calc SASA from Bio.PDB structure
693 #
694 # Usage
695 #
696 # result, sasa_classes = calcBioPDB(structure, ...)
697 #
698 # @remark Experimental, not thorougly tested yet
699 #
700 # @param bioPDBStructure A Bio.PDB structure
701 # @param parameters A freesasa.Paremeters object
702 # @param classifier A freesasa.Classifier object
703 # @param options Options supported are 'hetatm', 'skip-unknown' and 'halt-at-unknown'
704 #
705 # @return A freesasa.Result object and a dictionary with classes
706 # defined by the classifier and associated areas
707 #
708 # @exception Exception if unknown atom is encountered and the option
709 # 'halt-at-unknown' is active. Passes on exceptions from
710 # calc(), classifyResults() and structureFromBioPDB().
711 def calcBioPDB(bioPDBStructure, parameters = Parameters(),
712  classifier = None, options = Structure.defaultOptions):
713  structure = structureFromBioPDB(bioPDBStructure, classifier, options)
714  result = calc(structure, parameters)
715  sasa_classes = classifyResults(result, structure, classifier)
716  return result, sasa_classes
717 
718 
def structureFromBioPDB
Create a freesasa structure from a Bio.PDB structure.
Definition: freesasa.pyx:662
def setProbeRadius(self, r)
Set probe radius.
Definition: freesasa.pyx:109
const double * freesasa_structure_coord_array(const freesasa_structure *structure)
Get array of coordinates.
Assigns class and radius to atom by residue and atom name.
Definition: freesasa.pyx:222
Stores parameter values to be used by calculation.
Definition: freesasa.pyx:59
def calc
Calculate SASA of Structure.
Definition: freesasa.pyx:581
def __init__
Constructor.
Definition: freesasa.pyx:322
def probeRadius(self)
Get probe radius.
Definition: freesasa.pyx:115
def getVerbosity()
Get global verbosity.
Definition: freesasa.pyx:646
def __dealloc__(self)
The destructor.
Definition: freesasa.pyx:518
def selectArea(commands, structure, result)
Sum SASA result over a selection of atoms.
Definition: freesasa.pyx:620
def atomArea(self, i)
SASA for a given atom.
Definition: freesasa.pyx:199
def __init__
Initializes Parameters object.
Definition: freesasa.pyx:68
def nSlices(self)
Get the number of slices per atom in Lee & Richards algorithm.
Definition: freesasa.pyx:139
def structureArray
Create array of structures from PDB file.
Definition: freesasa.pyx:551
def radius(self, residueName, atomName)
Radius of atom.
Definition: freesasa.pyx:277
const char * freesasa_structure_atom_res_number(const freesasa_structure *structure, int i)
Get residue number.
Represents a protein structure, including its atomic radii.
Definition: freesasa.pyx:293
def addAtom(self, atomName, residueName, residueNumber, chainLabel, x, y, z)
Add atom to structure.
Definition: freesasa.pyx:368
char freesasa_structure_atom_chain(const freesasa_structure *structure, int i)
Get chain label.
freesasa_verbosity freesasa_get_verbosity(void)
Get the current verbosity level.
freesasa_structure * freesasa_structure_new(void)
Allocate empty structure.
freesasa_result * freesasa_calc_structure(const freesasa_structure *structure, const freesasa_parameters *parameters)
Calculates SASA based on a given structure.
int freesasa_structure_n(const freesasa_structure *structure)
Get number of atoms.
def setRadiiWithClassifier(self, classifier)
Assign radii to atoms in structure using a classifier.
Definition: freesasa.pyx:385
def radius(self, i)
Radius of atom.
Definition: freesasa.pyx:422
def nAtoms(self)
Number of atoms.
Definition: freesasa.pyx:414
def __cinit__(self)
The constructor.
Definition: freesasa.pyx:170
def calcBioPDB
Calc SASA from Bio.PDB structure.
Definition: freesasa.pyx:712
def atomName(self, i)
Get atom name.
Definition: freesasa.pyx:433
def setRadii(self, radiusArray)
Set atomic radii from an array.
Definition: freesasa.pyx:399
def nThreads(self)
Get the number of threads to use in calculations.
Definition: freesasa.pyx:153
const double * freesasa_structure_radius(const freesasa_structure *structure)
Returns a pointer to an array of the radii of each atom.
def classifyResults
Break SASA result down into classes.
Definition: freesasa.pyx:599
def coord(self, i)
Get coordinates of given atom.
Definition: freesasa.pyx:472
int freesasa_structure_add_atom(freesasa_structure *structure, const char *atom_name, const char *residue_name, const char *residue_number, char chain_label, double x, double y, double z)
Add individual atom to structure using default behavior.
def nAtoms(self)
Number of atoms in the results.
Definition: freesasa.pyx:180
def residueNumber(self, i)
Get residue number for given atom.
Definition: freesasa.pyx:451
def algorithm(self)
Get algorithm.
Definition: freesasa.pyx:99
void freesasa_classifier_free(freesasa_classifier *classifier)
Frees a classifier object.
def setVerbosity(verbosity)
Set global verbosity.
Definition: freesasa.pyx:639
int freesasa_select_area(const char *command, char *name, double *area, const freesasa_structure *structure, const freesasa_result *result)
Get area of a selection.
int freesasa_set_verbosity(freesasa_verbosity v)
Set the global verbosity level.
def __cinit__
Constructor.
Definition: freesasa.pyx:235
freesasa_classifier * freesasa_classifier_from_file(FILE *file)
Generate a classifier from a config-file.
def __dealloc__(self)
The destructor.
Definition: freesasa.pyx:249
def nPoints(self)
Get number of test points in Shrake & Rupley algorithm.
Definition: freesasa.pyx:127
def setNPoints(self, n)
Set number of test points in Shrake & Rupley algorithm.
Definition: freesasa.pyx:121
const char * freesasa_structure_atom_res_name(const freesasa_structure *structure, int i)
Get residue name.
def setAlgorithm(self, alg)
Set algorithm.
Definition: freesasa.pyx:88
Stores results from SASA calculation.
Definition: freesasa.pyx:166
def residueName(self, i)
Get residue name of given atom.
Definition: freesasa.pyx:442
void freesasa_structure_set_radius(freesasa_structure *structure, const double *radii)
Override the radii set when the structure was initialized.
def setNThreads(self, n)
Set the number of threads to use in calculations.
Definition: freesasa.pyx:146
freesasa_structure ** freesasa_structure_array(FILE *pdb, int *n, const freesasa_classifier *classifier, int options)
Init array of structures from PDB.
def classify(self, residueName, atomName)
Class of atom.
Definition: freesasa.pyx:262
def totalArea(self)
Total SASA.
Definition: freesasa.pyx:190
def chainLabel(self, i)
Get chain label for given atom.
Definition: freesasa.pyx:460
freesasa_structure * freesasa_structure_from_pdb(FILE *pdb, const freesasa_classifier *classifier, int options)
Init structure with coordinates from pdb-file.
def setNSlices(self, n)
Set the number of slices per atom in Lee & Richards algorithm.
Definition: freesasa.pyx:133
const char * freesasa_structure_atom_name(const freesasa_structure *structure, int i)
Get atom name.
void freesasa_structure_free(freesasa_structure *structure)
Free structure.
def __dealloc__(self)
The destructor.
Definition: freesasa.pyx:174
void freesasa_result_free(freesasa_result *result)
Frees a freesasa_result object.