FreeSASA  1.1
Open source SASA calculations
View on GitHub
Python interface

If Python is enabled using

$ ./configure --enable-python-bindings

Cython is used to build Python bindings for FreeSASA, and make install will install them.

Below follow some illustrations of how to use the package, see the package documentation for details.

Basic calculations

Using defaults everywhere a simple calculation can be carried out as follows (assuming the PDB structure 1UBQ is available)

1 import freesasa
2 
3 structure = freesasa.Structure("1ubq.pdb")
4 result = freesasa.calc(structure)
5 area_classes = freesasa.classifyResults(result,structure)
6 
7 print "Total : %.2f A2" % result.totalArea()
8 for key in area_classes:
9  print key, ": %.2f A2" % area_classes[key]

Which would give the following output

Total : 4804.06 A2
Polar : 2504.22 A2
Apolar : 2299.84 A2

The following does a high precision L&R calculation

1 result = freesasa.calc(structure,
2  freesasa.Parameters({'algorithm' : freesasa.LeeRichards,
3  'n-slices' : 100}))

Using the results from a calculation we can also integrate SASA over a selection of atoms, using a subset of the Pymol selection syntax (see Selection syntax):

1 selections = freesasa.selectArea(('alanine, resn ala','r1_10, resi 1-10'),
2  structure, result)
3 for key in selections:
4  print key, ": %.2f A2" % selections[key]

which gives the output

alanine : 120.08 A2
r1_10 : 634.31 A2

Customizing atom classification

This uses the NACCESS parameters (the file 'naccess.config' is available in the share/ directory of the repository).

1 classifier = freesasa.Classifier("naccess.config")
2 structure = freesasa.Structure("1ubq.pdb",classifier)
3 result = freesasa.calc(structure)
4 area_classes = freesasa.classifyResults(result,structure,classifier)

Classification can be customized also by extending the Classifier interface. The code below is an illustration of a classifier that classes Nitrogens separately, and assigns radii based on element only (and crudely).

1 import freesasa
2 import re
3 
4 class DerivedClassifier(Classifier):
5  def classify(self,residueName,atomName):
6  if re.match('\s*N',atomName):
7  return 'Nitrogen'
8  return 'Not-nitrogen'
9 
10  def radius(self,residueName,atomName):
11  if re.match('\s*N',atomName): # Nitrogen
12  return 1.6
13  if re.match('\s*C',atomName): # Carbon
14  return 1.7
15  if re.match('\s*O',atomName): # Oxygen
16  return 1.4
17  if re.match('\s*S',atomName): # Sulfur
18  return 1.8
19  return 0; # everything else (Hydrogen, etc)
20 
21 classifier = DerivedClassifier()
22 
23 # use the DerivedClassifier to calculate atom radii
24 structure = freesasa.Structure("1ubq.pdb",classifier)
25 result = freesasa.calc(structure)
26 
27 # use the DerivedClassifier to classify atoms
28 area_classes = freesasa.classifyResults(result,structure,classifier)

Of course, this example is somewhat contrived, if we only want the integrated area of Nitrogen atoms, the simpler choice would be

1 selection = freesasa.selectArea('nitrogen, symbol n', structure, result)

However, extending freesasa.Classifier, as illustrated above, allows classification to arbitrary complexity and also lets us redefine the radii used in the calculation.

Bio.PDB

FreeSASA can also calculate the SASA of a Bio.PDB structure

1 from Bio.PDB import PDBParser
2 parser = PDBParser()
3 structure = parser.get_structure("Ubiquitin","1ubq.pdb")
4 result, sasa_classes = freesasa.calcBioPDB(structure)

If one needs more control over the analysis the structure can be converted to a freesasa.Structure using freesasa.structureFromBioPDB() and the calculation can be performed the normal way using this structure.